MUSE Pipeline Reference Manual  0.18.5
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-2014 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 /*----------------------------------------------------------------------------*/
196 /*----------------------------------------------------------------------------*/
197 cpl_error_code
198 muse_wcs_locate_sources(muse_pixtable *aPixtable, float aSigma,
199  muse_wcs_centroid_type aCentroid,
200  muse_wcs_object *aWCSObj)
201 {
202  cpl_ensure_code(aPixtable && aPixtable->header && aWCSObj,
203  CPL_ERROR_NULL_INPUT);
204  cpl_ensure_code(aSigma >0., CPL_ERROR_ILLEGAL_INPUT);
205  switch (aCentroid) {
207  cpl_msg_info(__func__, "Gaussian profile fits for object centroiding");
208  break;
210  cpl_msg_info(__func__, "Moffat profile fits for object centroiding");
211  break;
213  cpl_msg_info(__func__, "Simple square box object centroiding");
214  break;
215  default:
216  return cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
217  }
218  if (getenv("MUSE_DEBUG_WCS") && atoi(getenv("MUSE_DEBUG_WCS")) > 2) {
219  const char *fn = "wcs__pixtable.fits";
220  cpl_msg_info(__func__, "Saving pixel table as \"%s\"", fn);
221  muse_pixtable_save(aPixtable, fn);
222  }
223 
224  /* check that DAR correction was carried out in some way */
225  cpl_boolean darcorrected = CPL_FALSE;
226  if (cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_DAR_NAME) &&
227  cpl_propertylist_get_double(aPixtable->header, MUSE_HDR_PT_DAR_NAME) > 0.) {
228  darcorrected = CPL_TRUE;
229  } else if (cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_DAR_CORR)) {
230  darcorrected = CPL_TRUE;
231  }
232  if (!darcorrected) {
233  const char *message = "Correction for differential atmospheric refraction "
234  "was not applied! Deriving astrometric correction "
235  "from such data is unlikely to give good results!";
236  cpl_msg_warning(__func__, "%s", message);
237  cpl_error_set_message(__func__, CPL_ERROR_UNSUPPORTED_MODE, "%s", message);
238  }
239 
240  muse_resampling_params *params =
242  params->pfx = 1.; /* large pixfrac to be sure to cover most gaps */
243  params->pfy = 1.;
244  params->pfl = 1.;
245  params->dlambda = 50.; /* we can integrate over lots of Angstrom here */
246  params->crtype = MUSE_RESAMPLING_CRSTATS_MEDIAN; /* median for very clean cube */
247  params->crsigma = 15.; /* very moderate CR rejection */
248  muse_datacube *cube = muse_resampling_cube(aPixtable, params, NULL);
250  /* reset cosmic ray statuses in aPixtable, since the CR rejection *
251  * done here might not be appropriate for the final datacube */
252  muse_pixtable_reset_dq(aPixtable, EURO3D_COSMICRAY);
253  if (!cube) {
254  return cpl_error_set_message(__func__, cpl_error_get_code(),
255  "Failure while creating cube!");
256  }
257  aWCSObj->cube = cube;
258  if (getenv("MUSE_DEBUG_WCS") && atoi(getenv("MUSE_DEBUG_WCS")) >= 2) {
259  const char *fn = "wcs__cube.fits";
260  cpl_msg_info(__func__, "Saving cube as \"%s\"", fn);
261  muse_datacube_save(cube, fn);
262  }
263 
264  /* detect objects in the cube, using image planes around the central one */
265  int iplane, irefplane = cpl_imagelist_get_size(cube->data) / 2;
266  /* medianing three images removes all cosmic rays but continuum objects stay *
267  * (need to use muse_images and the muse_combine function because *
268  * cpl_imagelist_collapse_median_create() disregards all bad pixels) */
270  unsigned int ilist = 0;
271  for (iplane = irefplane - 1; iplane <= irefplane + 1; iplane++) {
272  muse_image *image = muse_image_new();
273  image->data = cpl_image_duplicate(cpl_imagelist_get(cube->data, iplane));
274  image->dq = cpl_image_duplicate(cpl_imagelist_get(cube->dq, iplane));
275  image->stat = cpl_image_duplicate(cpl_imagelist_get(cube->stat, iplane));
276  muse_imagelist_set(list, image, ilist++);
277  } /* for iplane (planes around ref. wavelength) */
278  muse_image *median = muse_combine_median_create(list);
279  muse_imagelist_delete(list);
280  if (!median) {
281  return cpl_error_set_message(__func__, cpl_error_get_code(),
282  "Failure while creating detection image!");
283  }
284  if (getenv("MUSE_DEBUG_WCS") && atoi(getenv("MUSE_DEBUG_WCS")) >= 2) {
285  const char *fn = "wcs__image_median.fits";
286  cpl_msg_info(__func__, "Saving median detection image as \"%s\"", fn);
287  median->header = cpl_propertylist_new();
288  cpl_propertylist_append_string(median->header, "BUNIT",
289  cpl_propertylist_get_string(cube->header,
290  "BUNIT"));
291  muse_image_save(median, fn);
292  }
293  /* reject data in image to signify bad pixels to CPL routines */
295 
296  cpl_apertures *apertures = cpl_apertures_extract_sigma(median->data, aSigma);
297  int ndet = apertures ? cpl_apertures_get_size(apertures) : 0;
298  if (ndet < 1) {
299  muse_image_delete(median);
300  return cpl_error_set_message(__func__, CPL_ERROR_DATA_NOT_FOUND,
301  "No sources found for astrometric calibration "
302  "down to %.1f sigma limit on plane %d", aSigma,
303  irefplane + 1);
304  }
305  cpl_msg_debug(__func__, "The %f sigma threshold was used to find %d sources "
306  "in cube plane %d", aSigma, ndet, irefplane + 1);
307  /* save pixel value and sky coordinates of center of the data */
308  aWCSObj->xcenter = cpl_image_get_size_x(median->data) / 2.,
309  aWCSObj->ycenter = cpl_image_get_size_y(median->data) / 2.;
310  aWCSObj->ra = muse_pfits_get_ra(aPixtable->header);
311  aWCSObj->dec = muse_pfits_get_dec(aPixtable->header);
312 #if 0
313  cpl_msg_debug(__func__, "image size: %d x %d --> center %f, %f",
314  (int)cpl_image_get_size_x(median->data),
315  (int)cpl_image_get_size_y(median->data),
316  aWCSObj->xcenter, aWCSObj->ycenter);
317 #endif
318 
319  /* need error image not variance! */
320  cpl_image_power(median->stat, 0.5);
321  int nx = cpl_image_get_size_x(median->data),
322  ny = cpl_image_get_size_y(median->data);
323  /* parameters for the fits */
324  const double bgguess = cpl_image_get_median(median->data),
325  beta = 2.5, /* moffat beta */
326  moffat_alpha_to_fwhm = 1 / (2 * sqrt(pow(2, 1/beta) - 1));
327 #if 0 /* unused idea to further robustify FWHM estimate */
328  double fwhmguess = (muse_pfits_get_fwhm_start(aPixtable->header)
329  + muse_pfits_get_fwhm_end(aPixtable->header)) / 2.;
330  if (fwhmguess < FLT_EPSILON) { /* in case FWHM headers are not there */
331  fwhmguess = 4.; /* [pix] assume median seeing in WFM */
332  } else { /* use nominal MUSE spaxel scale to get FWHM guess [pix] */
333  if (muse_pfits_get_mode(aPixtable->header) < MUSE_MODE_NFM_AO_N) {
334  fwhmguess /= (kMuseSpaxelSizeX_WFM + kMuseSpaxelSizeX_WFM) / 2;
335  } else {
336  fwhmguess /= (kMuseSpaxelSizeX_NFM + kMuseSpaxelSizeX_NFM) / 2;
337  }
338  }
339 #endif
340 
341  cpl_table *detections = muse_cpltable_new(muse_wcs_detections_def, ndet);
342  cpl_table_unselect_all(detections);
343  int idet;
344  for (idet = 0; idet < ndet; idet++) {
345  double xc = cpl_apertures_get_centroid_x(apertures, idet+1),
346  yc = cpl_apertures_get_centroid_y(apertures, idet+1),
347  fluxaper = cpl_apertures_get_flux(apertures, idet+1);
348  double xwaper, ywaper;
349  cpl_image_get_fwhm(median->data, lround(xc), lround(yc), &xwaper, &ywaper);
350  /* Remove the detection if the FWHM could not be computed in both *
351  * axes, since this likely points to an object on the border or with *
352  * some other problem which should not be relied on for centroiding. */
353  if (xwaper < 0 || ywaper < 0) {
354 #if 1 /* XXX */
355  cpl_msg_debug(__func__, "FWHM computation failed at %f,%f, result was "
356  "%.3f,%.3f", xc, yc, xwaper, ywaper);
357 #endif
358  cpl_table_select_row(detections, idet);
359  continue;
360  }
361 
362  /* set the halfsize of the window to measure the centroids *
363  * with respect to the detection; +/-5 pix should suffice */
364  int x1 = floor(xc) - 5,
365  x2 = ceil(xc) + 5,
366  y1 = floor(yc) - 5,
367  y2 = ceil(yc) + 5;
368  /* force window to be inside the image */
369  if (x1 < 1) x1 = 1;
370  if (y1 < 1) y1 = 1;
371  if (x2 > nx) x2 = nx;
372  if (y2 > ny) y2 = ny;
373 
374  double xcen, ycen, xerr = 0., yerr = 0., xw, yw, flux;
375  cpl_errorstate state = cpl_errorstate_get();
376  switch (aCentroid) {
378  cpl_array *pgauss = cpl_array_new(7, CPL_TYPE_DOUBLE),
379  *pgerr = cpl_array_new(7, CPL_TYPE_DOUBLE),
380  *pgfit = cpl_array_new(7, CPL_TYPE_INT);
381  /* first, only fit the centroid */
382  cpl_array_set(pgfit, 3, 1); /* xc */
383  cpl_array_set(pgfit, 4, 1); /* yc */
384  cpl_array_set(pgfit, 0, 0);
385  cpl_array_set(pgfit, 1, 0);
386  cpl_array_set(pgfit, 2, 0);
387  cpl_array_set(pgfit, 5, 0);
388  cpl_array_set(pgfit, 6, 0);
389  cpl_array_set(pgauss, 3, xc);
390  cpl_array_set(pgauss, 4, yc);
391  cpl_array_set(pgauss, 0, bgguess);
392  cpl_array_set(pgauss, 1, fluxaper);
393  cpl_array_set(pgauss, 2, 0.); /* rho */
394  cpl_array_set(pgauss, 5, xwaper * CPL_MATH_SIG_FWHM);
395  cpl_array_set(pgauss, 6, ywaper * CPL_MATH_SIG_FWHM);
396  cpl_fit_image_gaussian(median->data, median->stat, lround(xc), lround(yc),
397  x2-x1+1, y2-y1+1, pgauss, pgerr, pgfit, NULL, NULL,
398  NULL, NULL, NULL, NULL, NULL);
399  xcen = cpl_array_get(pgauss, 3, NULL);
400  ycen = cpl_array_get(pgauss, 4, NULL);
401  xerr = cpl_array_get(pgerr, 3, NULL);
402  yerr = cpl_array_get(pgerr, 4, NULL);
403  /* second, keep center fixed, only fit flux and width */
404  cpl_array_set(pgfit, 3, 0); /* xc */
405  cpl_array_set(pgfit, 4, 0); /* yc */
406  cpl_array_set(pgfit, 1, 1); /* flux */
407  cpl_array_set(pgfit, 5, 1); /* sigma_x */
408  cpl_array_set(pgfit, 6, 1); /* sigma_y */
409  cpl_fit_image_gaussian(median->data, median->stat, lround(xc), lround(yc),
410  x2-x1+1, y2-y1+1, pgauss, pgerr, pgfit, NULL, NULL,
411  NULL, NULL, NULL, NULL, NULL);
412  xw = cpl_array_get(pgauss, 5, NULL) * CPL_MATH_FWHM_SIG;
413  yw = cpl_array_get(pgauss, 6, NULL) * CPL_MATH_FWHM_SIG;
414  flux = cpl_array_get(pgauss, 1, NULL);
415  cpl_array_delete(pgauss);
416  cpl_array_delete(pgerr);
417  cpl_array_delete(pgfit);
418  break;
419  }
421  cpl_size nmax = (x2-x1+1) * (y2-y1+1);
422  cpl_matrix *pos = cpl_matrix_new(nmax, 2);
423  cpl_vector *val = cpl_vector_new(nmax),
424  *err = cpl_vector_new(nmax);
425  int i, npix = 0;
426  for (i = x1; i <= x2; i++) {
427  int j;
428  for (j = y1; j <= y2; j++) {
429  int error;
430  double value = cpl_image_get(median->data, i, j, &error);
431  if (error != 0) {
432  continue;
433  }
434  cpl_matrix_set(pos, npix, 0, i);
435  cpl_matrix_set(pos, npix, 1, j);
436  cpl_vector_set(val, npix, value);
437  /* stat is already sigma! */
438  cpl_vector_set(err, npix, cpl_image_get(median->stat, i, j, &error));
439  npix++;
440  } /* for j */
441  } /* for i */
442  cpl_matrix_set_size(pos, npix, 2);
443  cpl_vector_set_size(val, npix);
444  cpl_vector_set_size(err, npix);
445  cpl_array *pmoffat = cpl_array_new(8, CPL_TYPE_DOUBLE),
446  *pmerror = cpl_array_new(8, CPL_TYPE_DOUBLE),
447  *pmfit = cpl_array_new(8, CPL_TYPE_INT);
448  /* first, only fit the centroid */
449  cpl_array_set(pmfit, 2, 1); /* xc */
450  cpl_array_set(pmfit, 3, 1); /* yc */
451  cpl_array_set(pmfit, 0, 0);
452  cpl_array_set(pmfit, 1, 0);
453  cpl_array_set(pmfit, 4, 0);
454  cpl_array_set(pmfit, 5, 0);
455  cpl_array_set(pmfit, 6, 0);
456  cpl_array_set(pmfit, 7, 0);
457  cpl_array_set(pmoffat, 2, xc);
458  cpl_array_set(pmoffat, 3, yc);
459  cpl_array_set(pmoffat, 0, bgguess);
460  cpl_array_set(pmoffat, 1, fluxaper);
461  cpl_array_set(pmoffat, 4, xwaper * moffat_alpha_to_fwhm);
462  cpl_array_set(pmoffat, 5, ywaper * moffat_alpha_to_fwhm);
463  cpl_array_set(pmoffat, 6, beta);
464  cpl_array_set(pmoffat, 7, 0.); /* rho */
465  muse_utils_fit_moffat_2d(pos, val, err, pmoffat, pmerror, pmfit, NULL, NULL);
466  xcen = cpl_array_get(pmoffat, 2, NULL);
467  ycen = cpl_array_get(pmoffat, 3, NULL);
468  xerr = cpl_array_get(pmerror, 2, NULL);
469  yerr = cpl_array_get(pmerror, 3, NULL);
470  /* second, keep center fixed, only fit flux and width */
471  cpl_array_set(pmfit, 2, 0); /* xc */
472  cpl_array_set(pmfit, 3, 0); /* yc */
473  cpl_array_set(pmfit, 1, 1); /* flux */
474  cpl_array_set(pmfit, 4, 1); /* xhwhm */
475  cpl_array_set(pmfit, 5, 1); /* yhwhm */
476  muse_utils_fit_moffat_2d(pos, val, err, pmoffat, pmerror, pmfit, NULL, NULL);
477  xw = cpl_array_get(pmoffat, 4, NULL) / moffat_alpha_to_fwhm;
478  yw = cpl_array_get(pmoffat, 5, NULL) / moffat_alpha_to_fwhm;
479  flux = cpl_array_get(pmoffat, 1, NULL);
480  cpl_array_delete(pmoffat);
481  cpl_array_delete(pmerror);
482  cpl_array_delete(pmfit);
483  cpl_matrix_delete(pos);
484  cpl_vector_delete(val);
485  cpl_vector_delete(err);
486  break;
487  }
488  default: /* MUSE_WCS_CENTROID_BOX is left */
489  muse_utils_image_get_centroid_window(median->data, x1, y1, x2, y2,
490  &xcen, &ycen,
492 #if 0
493  printf("%d apertures %f %f boxes %f %f deltas %f %f\n", idet+1, xc, yc,
494  xcen, ycen, xcen - xc, ycen - yc);
495  fflush(stdout);
496 #endif
497  /* set error to 0.15 pix which is the typical *
498  * stdev compared to Gaussian fits */
499  xerr = 0.15;
500  yerr = 0.15;
501  /* compute FWHM again with the revised central position */
502  cpl_image_get_fwhm(median->data, lround(xcen), lround(ycen), &xw, &yw);
503  flux = fluxaper; /* take the aperture flux */
504  }
505 
506  /* mandatory columns: position */
507  cpl_table_set(detections, "XPOS", idet, xcen);
508  cpl_table_set(detections, "YPOS", idet, ycen);
509  /* and the error on the position */
510  cpl_table_set(detections, "XERR", idet, xerr);
511  cpl_table_set(detections, "YERR", idet, yerr);
512  /* extra columns: FWHM... */
513  cpl_table_set(detections, "XFWHM", idet, xw);
514  cpl_table_set(detections, "YFWHM", idet, yw);
515  /* ... and flux */
516  cpl_table_set(detections, "FLUX", idet, flux);
517 
518  if (cpl_errorstate_is_equal(state) && xw > 0 && yw > 0 &&
519  xcen >= 1 && xcen <= nx && ycen >= 1 && ycen <= ny) {
520  continue;
521  }
522  /* some error occurred, so mark this entry for removal */
523  cpl_table_select_row(detections, idet);
524  } /* for idet (aperture index) */
525  cpl_table_erase_selected(detections);
526  cpl_apertures_delete(apertures);
527  muse_image_delete(median);
528  ndet = cpl_table_get_nrow(detections);
529  cpl_msg_debug(__func__, "%d valid sources were recorded in the detections "
530  "table", ndet);
531 
532  if (getenv("MUSE_DEBUG_WCS") && atoi(getenv("MUSE_DEBUG_WCS")) >= 2) {
533  const char *fn = "wcs__detections.fits";
534  cpl_msg_info(__func__, "Saving table of detections as \"%s\"", fn);
535  cpl_table_save(detections, NULL, NULL, fn, CPL_IO_CREATE);
536  }
537 
538  aWCSObj->detected = detections;
539  return CPL_ERROR_NONE;
540 } /* muse_wcs_locate_sources() */
541 
542 /*----------------------------------------------------------------------------*/
580 /*----------------------------------------------------------------------------*/
581 cpl_error_code
582 muse_wcs_solve(muse_wcs_object *aWCSObj, cpl_table *aReference,
583  float aRadius, float aFAccuracy, int aIter, float aThresh)
584 {
585  cpl_ensure_code(aWCSObj, CPL_ERROR_NULL_INPUT);
586  cpl_table *detected = aWCSObj->detected;
587  int ndet = cpl_table_get_nrow(detected),
588  nref = cpl_table_get_nrow(aReference);
589  cpl_ensure_code(ndet > 0 && nref > 0, CPL_ERROR_ILLEGAL_INPUT);
590  cpl_ensure_code(muse_cpltable_check(detected, muse_wcs_detections_def)
591  == CPL_ERROR_NONE &&
592  muse_cpltable_check(aReference, muse_wcs_reference_def)
593  == CPL_ERROR_NONE,
594  CPL_ERROR_BAD_FILE_FORMAT);
595  cpl_ensure_code(aRadius > 0.&& aFAccuracy > 0., CPL_ERROR_ILLEGAL_INPUT);
596 
597  /* sort tables by the source brightness */
598  cpl_propertylist *order = cpl_propertylist_new();
599  cpl_propertylist_append_bool(order, "FLUX", TRUE);
600  cpl_table_sort(detected, order);
601  cpl_propertylist_erase(order, "FLUX");
602  cpl_propertylist_append_bool(order, "mag", FALSE);
603  cpl_table_sort(aReference, order);
604  cpl_propertylist_delete(order);
605  if (getenv("MUSE_DEBUG_WCS") && atoi(getenv("MUSE_DEBUG_WCS")) >= 2) {
606  FILE *fp = fopen("wcs__detections.ascii", "w");
607  fprintf(fp, "%s: table of detected sources (sorted by flux):\n", __func__);
608  cpl_table_dump(detected, 0, 1000, fp);
609  fclose(fp);
610  fp = fopen("wcs__references.ascii", "w");
611  fprintf(fp, "%s: table of reference objects (sorted by flux):\n", __func__);
612  cpl_table_dump(aReference, 0, 1000, fp);
613  fclose(fp);
614  }
615 
616  /* construct a basic input WCS */
617  cpl_propertylist *wcsin = muse_wcs_create_default();
618  cpl_propertylist_append_double(wcsin, "CRVAL1", aWCSObj->ra);
619  cpl_propertylist_append_double(wcsin, "CRVAL2", aWCSObj->dec);
620  cpl_propertylist_update_double(wcsin, "CRPIX1", aWCSObj->crpix1);
621  cpl_propertylist_update_double(wcsin, "CRPIX2", aWCSObj->crpix2);
622  /* add NAXIS to fool cpl_wcs_new_from_propertylist() */
623  cpl_propertylist_append_int(wcsin, "NAXIS", 2);
624  cpl_propertylist_append_int(wcsin, "NAXIS1", aWCSObj->xcenter * 2.);
625  cpl_propertylist_append_int(wcsin, "NAXIS2", aWCSObj->ycenter * 2.);
626 
627  /* convert input tables into matrices for the pattern-matching function */
628  cpl_matrix *data = cpl_matrix_new(2, ndet),
629  *patt = cpl_matrix_new(2, nref);
630  int i;
631  for (i = 0; i < ndet; i++) {
632  cpl_matrix_set(data, 0, i, cpl_table_get(detected, "XPOS", i, NULL));
633  cpl_matrix_set(data, 1, i, cpl_table_get(detected, "YPOS", i, NULL));
634  } /* for i (all detections) */
635 
636  /* use the basic WCS to transform input reference positions *
637  * to pixel positions, to take out deformations by gnomonic *
638  * projection before attempting pattern matching */
639  for (i = 0; i < nref; i++) {
640  double ra = cpl_table_get(aReference, "RA", i, NULL),
641  dec = cpl_table_get(aReference, "DEC", i, NULL),
642  x, y;
643  muse_wcs_pixel_from_celestial(wcsin, ra, dec, &x, &y);
644  cpl_matrix_set(patt, 0, i, x);
645  cpl_matrix_set(patt, 1, i, y);
646 #if 0
647  printf("conversion: %2d\t%.7f %.7f\t%6.2f %6.2f\n", i + 1, ra, dec, x, y);
648  fflush(stdout);
649 #endif
650  } /* for i (all reference points) */
651 #if 0
652  if (cpl_msg_get_level() == CPL_MSG_DEBUG) {
653  printf("%s: data matrix:\n", __func__);
654  cpl_matrix_dump(data, stdout);
655  printf("%s: patt matrix:\n", __func__);
656  cpl_matrix_dump(patt, stdout);
657  fflush(stdout);
658  }
659 #endif
660 
661  /* compute typical data error from the error columns, *
662  * to use this as input for the pattern matching */
663  double accuracy = sqrt(pow(cpl_table_get_column_mean(detected, "XERR"), 2) +
664  pow(cpl_table_get_column_mean(detected, "YERR"), 2));
665  /* start with somewhat worse accuracy to make *
666  * it more likely that a match is found at all */
667  accuracy *= aFAccuracy;
668  double radius = aRadius; /* start with pattern matching of 5 pix search radius */
669  int nid = INT_MAX; /* number of identified detections */
670  cpl_array *aid = NULL;
671  cpl_boolean dupes = CPL_FALSE;
672  double xscale, yscale;
673  do {
674  double scale, angle;
675  /* As recommended in the CPL docs, initially select the 20 first (brightest!) *
676  * detections and the 10 first (brightest) reference sources, so that *
677  * hopefully all reference sources are contained in the brightest detections. */
678 #define USE_DATA 20
679 #define USE_PATT 10
680  int ndata = ndet < USE_DATA ? ndet : USE_DATA,
681  npatt = nref < USE_PATT ? nref : USE_PATT;
682  cpl_array_delete(aid);
683  do {
684  cpl_msg_debug(__func__, "trying pattern matching with accuracy %g and radius %g",
685  accuracy, radius);
686  aid = cpl_ppm_match_points(data, ndata, accuracy,
687  patt, npatt, 0.01 /* [pix] almost exact */,
688  0.1 /* 10% */, radius, NULL, NULL,
689  &scale, &angle);
690  /* decrease accuracy in case pattern matching didn't succeed */
691  accuracy /= aid ? 1. : 1.5;
692  if (accuracy < FLT_EPSILON) {
693  break; /* doesn't make sense to go any lower */
694  }
695  } while (!aid); /* no matched points likely means to low accuracy */
696  cpl_errorstate state = cpl_errorstate_get();
697  nid = cpl_array_get_size(aid);
698  if (nid > 0) { /* subtract valid only if the array exists */
699  nid -= cpl_array_count_invalid(aid);
700  }
701 #if 0
702  printf("aid (valid=%d):\n", nid);
703  cpl_array_dump(aid, 0, cpl_array_get_size(aid), stdout);
704  fflush(stdout);
705 #endif
706  if (nid < 1) {
707  cpl_array_delete(aid);
708  cpl_matrix_delete(data);
709  cpl_matrix_delete(patt);
710  cpl_errorstate_set(state); /* swallow error about NULL cpl_array */
711  cpl_propertylist_delete(wcsin);
712  return cpl_error_set_message(__func__, CPL_ERROR_DATA_NOT_FOUND, "None of "
713  "the %d detections could be identified with "
714  "the %d reference positions with radius %.3f "
715  "pix and data accuracy %.3e pix",
716  ndet, nref, radius, accuracy);
717  }
718  /* the first-guess scale computed here only makes sense *
719  * if multiplied by the pixel scale in the first-guess WCS */
720  muse_wcs_get_scales(wcsin, &xscale, &yscale);
721  dupes = muse_cplarray_has_duplicate(aid);
722  cpl_msg_debug(__func__, "%d %sidentified points (scale = %g, angle = %g; "
723  "used radius = %.3f pix, data accuracy = %.3e pix)", nid,
724  dupes ? "(non-unique!) " : "unique ",
725  scale*1800.*(xscale+yscale), angle, radius, accuracy);
726  radius /= 1.5; /* next loop with smaller radius, if necessary */
727  } while (dupes);
728  cpl_matrix_delete(data);
729  cpl_matrix_delete(patt);
730 
731  /* create matrices again for cpl_wcs_platesol(), this *
732  * time with detected pixel positions and sky positions, *
733  * but only for the identified detections */
734  cpl_matrix *mpx = cpl_matrix_new(nid, 2),
735  *msky = cpl_matrix_new(nid, 2);
736  int iid = 0; /* index of identified object in matrices */
737  for (i = 0; i < cpl_array_get_size(aid); i++) { /* index in reference table */
738  if (!cpl_array_is_valid(aid, i)) {
739  continue;
740  }
741  int idata = cpl_array_get_int(aid, i, NULL); /* index in detections table */
742  cpl_matrix_set(mpx, iid, 0, cpl_table_get(detected, "XPOS", idata, NULL));
743  cpl_matrix_set(mpx, iid, 1, cpl_table_get(detected, "YPOS", idata, NULL));
744  cpl_matrix_set(msky, iid, 0, cpl_table_get(aReference, "RA", i, NULL));
745  cpl_matrix_set(msky, iid, 1, cpl_table_get(aReference, "DEC", i, NULL));
746 #if 0
747  printf("matrices: %2d\t%.7f %.7f\t%6.2f %6.2f\n", iid + 1,
748  cpl_table_get(aReference, "RA", i, NULL),
749  cpl_table_get(aReference, "DEC", i, NULL),
750  cpl_table_get(detected, "XPOS", idata, NULL),
751  cpl_table_get(detected, "YPOS", idata, NULL));
752 #endif
753  iid++;
754  }
755 #if 0
756  printf("mpx:\n");
757  cpl_matrix_dump(mpx, stdout);
758  printf("msky:\n");
759  cpl_matrix_dump(msky, stdout);
760  fflush(stdout);
761 #endif
762  cpl_array_delete(aid);
763 
764  /* compute the solution */
765  cpl_propertylist *wcsout = NULL;
766  cpl_error_code rc = cpl_wcs_platesol(wcsin, msky, mpx, aIter, aThresh,
767  CPL_WCS_PLATESOL_6, CPL_WCS_MV_CRVAL,
768  &wcsout);
769  if (aWCSObj->cube) {
770  cpl_propertylist_copy_property_regexp(wcsout, aWCSObj->cube->header,
771  CPL_WCS_REGEXP, 1);
772  } /* if cube */
773  cpl_matrix_delete(mpx);
774  cpl_matrix_delete(msky);
775  cpl_propertylist_delete(wcsin);
776  if (rc != CPL_ERROR_NONE) {
777  cpl_msg_warning(__func__, "Computing the WCS solution returned an error "
778  "(%d): %s", rc, cpl_error_get_message());
779  }
780 
781  /* Compute the rotation angle and scales */
782  double cd11 = cpl_propertylist_get_double(wcsout, "CD1_1"),
783  cd22 = cpl_propertylist_get_double(wcsout, "CD2_2"),
784  cd12 = cpl_propertylist_get_double(wcsout, "CD1_2"),
785  cd21 = cpl_propertylist_get_double(wcsout, "CD2_1"),
786  det = cd11 * cd22 - cd12 * cd21;
787  if (det < 0.) {
788  cd12 *= -1;
789  cd11 *= -1;
790  }
791  /* the values we want, by default for non-rotation */
792  double xang, yang;
793  muse_wcs_get_angles(wcsout, &xang, &yang);
794  muse_wcs_get_scales(wcsout, &xscale, &yscale);
795  xscale *= 3600.; /* scales in arc seconds */
796  yscale *= 3600.;
797  cpl_msg_info(__func__, "astrometric calibration results: scales %f/%f "
798  "arcsec/spaxel, rotation %g/%g deg", xscale, yscale, xang, yang);
799 
800 #if 0
801  printf("%s: output propertylist (rc = %d):\n", __func__, rc);
802  fflush(stdout);
803  cpl_propertylist_save(wcsout, "astrometry_wcsout.fits", CPL_IO_CREATE);
804  system("fold -w 80 astrometry_wcsout.fits | grep -v \"^ \"");
805  remove("astrometry_wcsout.fits");
806 #endif
807 
808  /* number of stars input to the astrometric fit as QC */
809  cpl_propertylist_update_int(wcsout, QC_ASTROMETRY_NSTARS, nid);
810  /* scales for QC in arcsec */
811  cpl_propertylist_update_float(wcsout, QC_ASTROMETRY_SCX, xscale);
812  cpl_propertylist_update_float(wcsout, QC_ASTROMETRY_SCY, yscale);
813  /* angles in degrees */
814  cpl_propertylist_update_float(wcsout, QC_ASTROMETRY_ANGX, xang);
815  cpl_propertylist_update_float(wcsout, QC_ASTROMETRY_ANGY, yang);
816  /* copy the "systematic error" as residuals for QC */
817  double medresx = cpl_propertylist_get_double(wcsout, "CSYER1") * 3600.,
818  medresy = cpl_propertylist_get_double(wcsout, "CSYER2") * 3600.;
819  cpl_propertylist_update_float(wcsout, QC_ASTROMETRY_RESX, medresx);
820  cpl_propertylist_update_float(wcsout, QC_ASTROMETRY_RESY, medresy);
821 
822  aWCSObj->wcs = wcsout;
823  return rc;
824 } /* muse_wcs_solve() */
825 
826 /*----------------------------------------------------------------------------*/
836 /*----------------------------------------------------------------------------*/
837 cpl_propertylist *
839 {
840  cpl_propertylist *wcs = cpl_propertylist_new();
841 
842  /* We only care about the spatial axes, pretend that we have a 2D image; *
843  * set WCSAXES to make fitsverify happy when it finds the other ^C headers. */
844  cpl_propertylist_append_int(wcs, "WCSAXES", 2);
845 
846  /* Check, if we are dealing with wcslib >= 4.23, which has the fix with the *
847  * floating point numbers. If an older version is found, use a small number *
848  * instead of zero for CRPIX, so that on loading CPL sees it as float value. */
849  double smallvalue = FLT_MIN;
850  const char *cpldesc = cpl_get_description(CPL_DESCRIPTION_DEFAULT); /* never fails */
851  /* search for the WCSLIB version string */
852  char *pcpldesc = strstr(cpldesc, "WCSLIB = ");
853  if (pcpldesc) {
854  pcpldesc += 8;
855  double wcslibversion = atof(pcpldesc);
856  if (wcslibversion >= 4.23) {
857  smallvalue = 0.;
858  } /* if newer wcslib */
859  cpl_msg_debug(__func__, "found wcslib %.2f, using smallvalue = %e",
860  wcslibversion, smallvalue);
861  } /* if wcslib string found */
862 
863  /* axis 1 */
864  /* CRPIX is the zero order correction to the astrometry, set zero here */
865  cpl_propertylist_append_double(wcs, "CRPIX1", smallvalue);
866  /* negative value, to signify that east is to the left */
867  cpl_propertylist_append_double(wcs, "CD1_1", -kMuseSpaxelSizeX_WFM / 3600);
868  cpl_propertylist_append_string(wcs, "CTYPE1", "RA---TAN");
869  cpl_propertylist_append_string(wcs, "CUNIT1", "deg");
870 
871  /* axis 2 */
872  cpl_propertylist_append_double(wcs, "CRPIX2", smallvalue);
873  cpl_propertylist_append_double(wcs, "CD2_2", kMuseSpaxelSizeY_WFM / 3600);
874  cpl_propertylist_append_string(wcs, "CTYPE2", "DEC--TAN");
875  cpl_propertylist_append_string(wcs, "CUNIT2", "deg");
876 
877  /* cross-terms are assumed to be not present (no rotation known!) */
878  cpl_propertylist_append_double(wcs, "CD1_2", 0.);
879  cpl_propertylist_append_double(wcs, "CD2_1", 0.);
880 
881  /* leave out CRVAL1/2, as we use this only to do the gnomonic projection */
882 
883  return wcs;
884 } /* muse_wcs_create_default() */
885 
886 /*----------------------------------------------------------------------------*/
901 /*----------------------------------------------------------------------------*/
902 cpl_propertylist *
903 muse_wcs_apply_cd(const cpl_propertylist *aExpHeader,
904  const cpl_propertylist *aCalHeader)
905 {
906  cpl_ensure(aCalHeader && aExpHeader, CPL_ERROR_NULL_INPUT, NULL);
907 
908  /* duplicate the input WCS because we want to adapt it to the current data */
909  cpl_propertylist *header = cpl_propertylist_duplicate(aCalHeader);
910 
911  /* Find field rotation (in radians) and create a CD matrix that reflects *
912  * the exposure position angle corrected by the astrometric solution. */
913  double pa = muse_astro_posangle(aExpHeader) * CPL_MATH_RAD_DEG,
914  xang, yang, xsc, ysc;
915  muse_wcs_get_scales(header, &xsc, &ysc);
916  muse_wcs_get_angles(header, &xang, &yang);
917  cpl_msg_debug(__func__, "WCS solution: scales %f / %f arcsec, angles %f / %f "
918  "deg", xsc * 3600., ysc * 3600., xang, yang);
919  /* create PCi_j matrix from the CDi_j in the header and invert it */
920  cpl_matrix *pc = cpl_matrix_new(2, 2);
921  cpl_matrix_set(pc, 0, 0, cpl_propertylist_get_double(header, "CD1_1") / xsc);
922  cpl_matrix_set(pc, 0, 1, cpl_propertylist_get_double(header, "CD1_2") / xsc);
923  cpl_matrix_set(pc, 1, 0, cpl_propertylist_get_double(header, "CD2_1") / ysc);
924  cpl_matrix_set(pc, 1, 1, cpl_propertylist_get_double(header, "CD2_2") / ysc);
925  cpl_matrix *pcinv = cpl_matrix_invert_create(pc);
926  cpl_matrix_delete(pc);
927  /* now create corrective CDi_j elements from the inverted PCi_j matrix */
928  double cd11cor, cd12cor, cd21cor, cd22cor;
929  if (pcinv) {
930  cd11cor = xsc * cpl_matrix_get(pcinv, 0, 0);
931  cd12cor = xsc * cpl_matrix_get(pcinv, 0, 1);
932  cd21cor = ysc * cpl_matrix_get(pcinv, 1, 0);
933  cd22cor = ysc * cpl_matrix_get(pcinv, 1, 1);
934  cpl_matrix_delete(pcinv);
935  } else {
936  cpl_msg_warning(__func__, "CD matrix of astrometric solution could not "
937  "be inverted! Assuming negligible zeropoint rotation.");
938  cd11cor = xsc * 1.;
939  cd12cor = xsc * 0.;
940  cd21cor = ysc * 0.;
941  cd22cor = ysc * 1.;
942  }
943  /* now we can finally compute the effective CDi_j of the exposure */
944  double cd11 = cos(pa) * cd11cor - sin(pa) * cd21cor,
945  cd12 = cos(pa) * cd12cor - sin(pa) * cd22cor,
946  cd21 = sin(pa) * cd11cor + cos(pa) * cd21cor,
947  cd22 = sin(pa) * cd12cor + cos(pa) * cd22cor;
948  cpl_propertylist_update_double(header, "CD1_1", cd11),
949  cpl_propertylist_update_double(header, "CD1_2", cd12),
950  cpl_propertylist_update_double(header, "CD2_1", cd21);
951  cpl_propertylist_update_double(header, "CD2_2", cd22),
952  muse_wcs_get_scales(header, &xsc, &ysc);
953  muse_wcs_get_angles(header, &xang, &yang);
954  cpl_msg_debug(__func__, "Updated CD matrix (%f deg field rotation): "
955  "%e %e %e %e (scales %f / %f arcsec, angles %f / %f deg)",
956  pa * CPL_MATH_DEG_RAD, cd11, cd12, cd21, cd22,
957  xsc * 3600., ysc * 3600., xang, yang);
958  return header;
959 } /* muse_wcs_apply_cd() */
960 
961 /*----------------------------------------------------------------------------*/
990 /*----------------------------------------------------------------------------*/
991 cpl_error_code
992 muse_wcs_project_tan(muse_pixtable *aPixtable, const cpl_propertylist *aWCS)
993 {
994  cpl_size nrow = muse_pixtable_get_nrow(aPixtable);
995  /* nrow == 0 implies a NULL or broken pixel table */
996  cpl_ensure_code(nrow > 0 && aPixtable->header && aWCS, CPL_ERROR_NULL_INPUT);
997  cpl_ensure_code(muse_pixtable_wcs_check(aPixtable) == MUSE_PIXTABLE_WCS_PIXEL,
998  CPL_ERROR_INCOMPATIBLE_INPUT);
999  /* ensure that the input WCS is for gnomonic projection */
1000  const char *type1 = cpl_propertylist_get_string(aWCS, "CTYPE1"),
1001  *type2 = cpl_propertylist_get_string(aWCS, "CTYPE2");
1002  cpl_ensure_code(type1 && type2 && !strncmp(type1, "RA---TAN", 9) &&
1003  !strncmp(type2, "DEC--TAN", 9),
1004  CPL_ERROR_UNSUPPORTED_MODE);
1005 
1006  /* clean main WCS keys from input pixel table, in *
1007  * case there are keys we do not overwrite below */
1008  cpl_propertylist_erase_regexp(aPixtable->header, MUSE_WCS_KEYS, 0);
1009 
1010  /* apply the exposure rotation on top of the zeropoint *
1011  * rotation from the astrometric calibration */
1012  cpl_propertylist *header = muse_wcs_apply_cd(aPixtable->header, aWCS);
1013  /* don't want CRVAL or LON/LATPOLE in the output, *
1014  * because we create native coords here */
1015  cpl_propertylist_erase_regexp(header, "^CRVAL[0-9]+$|^L[OA][NT]POLE$", 0);
1016  double cd11 = cpl_propertylist_get_double(header, "CD1_1"),
1017  cd12 = cpl_propertylist_get_double(header, "CD1_2"),
1018  cd21 = cpl_propertylist_get_double(header, "CD2_1"),
1019  cd22 = cpl_propertylist_get_double(header, "CD2_2");
1020 
1021  /* Compute reference pixel from center of the data plus the zero *
1022  * order correction of the corrective WCS (CRPIX1/2 of aWCS); use *
1023  * the spatial extents before the DAR correction, if possible. */
1024  cpl_errorstate prestate = cpl_errorstate_get();
1025  double xlo = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_PREDAR_XLO),
1026  xhi = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_PREDAR_XHI),
1027  ylo = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_PREDAR_YLO),
1028  yhi = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_PREDAR_YHI);
1029  if (!cpl_errorstate_is_equal(prestate)) {
1030  cpl_errorstate_set(prestate);
1031  /* try normal pixel table headers now */
1032  xlo = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_XLO);
1033  xhi = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_XHI);
1034  ylo = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_YLO);
1035  yhi = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_YHI);
1036  }
1037  double wcspix1 = cpl_propertylist_get_double(header, "CRPIX1"),
1038  wcspix2 = cpl_propertylist_get_double(header, "CRPIX2"),
1039  crpix1 = (xhi + xlo) / 2. + wcspix1,
1040  crpix2 = (yhi + ylo) / 2. + wcspix2;
1041  cpl_propertylist_update_double(header, "CRPIX1", crpix1),
1042  cpl_propertylist_update_double(header, "CRPIX2", crpix2),
1043  cpl_msg_debug(__func__, "Using reference pixel %f/%f (limits in pixel table "
1044  "%f..%f/%f..%f, WCS correction %f,%f)", crpix1, crpix2,
1045  xlo, xhi, ylo, yhi, wcspix1, wcspix2);
1046 
1047  /* delete the units of the x/y columns while we are working on them */
1048  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_XPOS, "");
1049  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_YPOS, "");
1050 
1051  /* replace coordinate values in the pixel table by the computed ones, *
1052  * carry out the _partial_ WCS coordinate transform for each pixel */
1053  float *xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
1054  *ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS);
1055  cpl_size n;
1056  #pragma omp parallel for default(none) /* as req. by Ralf */ \
1057  shared(cd11, cd12, cd21, cd22, crpix1, crpix2, nrow, xpos, ypos)
1058  for (n = 0; n < nrow; n++) {
1059  /* conversion from pixel coordinates to projection plane coordinates; *
1060  * x,y as in Calabretta & Greisen 2002 A&A 395, 1077 (Paper II), Fig. 1 */
1061  double x = cd11 * (xpos[n] - crpix1) + cd12 * (ypos[n] - crpix2),
1062  y = cd22 * (ypos[n] - crpix2) + cd21 * (xpos[n] - crpix1);
1063 
1064  /* conversion from projection plane to native spherical coordinates: *
1065  * phi,theta as in Calabretta & Greisen 2002 (Paper II), Fig. 1; *
1066  * these formulae are for the gnomomic (TAN) case, Sect. 5.1/5.1.3, *
1067  * i.e. Eq. (14), (15), and (55) *
1068  * As we only further use these in other parts of the pipeline, the *
1069  * values are not converted to degrees but stay in radians. */
1070  double phi = atan2(x, -y),
1071  theta = atan(CPL_MATH_DEG_RAD / sqrt(x*x + y*y));
1072  if (phi < 0) { /* phi should be between 0 and 2pi, to let tests pass */
1073  phi += CPL_MATH_2PI;
1074  }
1075 
1076  /* conversion from native spherical to celestial spherical coordinates *
1077  * is done later when combining exposures, see muse_xcombine_tables() */
1078  xpos[n] = phi;
1079  ypos[n] = theta - CPL_MATH_PI_2; /* subtract pi/2 for better accuracy */
1080  } /* for n (table/matrix rows) */
1081 
1082  /* Here we produce units in radians; use "rad" unit string to signify *
1083  * that these are native spherical coordinates but haven't gotten full *
1084  * WCS treatment to alpha,delta yet */
1085  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_XPOS, "rad");
1086  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_YPOS, "rad");
1087  muse_pixtable_compute_limits(aPixtable);
1088  /* copy spatial WCS to the pixel table */
1089  cpl_propertylist_copy_property_regexp(aPixtable->header, header,
1090  MUSE_WCS_KEYS, 0);
1091  cpl_propertylist_delete(header);
1092 
1093  /* add the status header */
1094  cpl_propertylist_update_string(aPixtable->header, MUSE_HDR_PT_WCS,
1096  cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_WCS,
1097  MUSE_HDR_PT_WCS_COMMENT_PROJ);
1098  return CPL_ERROR_NONE;
1099 } /* muse_wcs_project_tan() */
1100 
1101 /*----------------------------------------------------------------------------*/
1130 /*----------------------------------------------------------------------------*/
1131 cpl_error_code
1132 muse_wcs_position_celestial(muse_pixtable *aPixtable, double aRA, double aDEC)
1133 {
1134  cpl_size nrow = muse_pixtable_get_nrow(aPixtable);
1135  /* nrow == 0 implies a NULL or broken pixel table */
1136  cpl_ensure_code(nrow > 0 && aPixtable->header, CPL_ERROR_NULL_INPUT);
1137  cpl_ensure_code(muse_pixtable_wcs_check(aPixtable) == MUSE_PIXTABLE_WCS_NATSPH,
1138  CPL_ERROR_INCOMPATIBLE_INPUT);
1139  /* ensure that the input WCS is for gnomonic projection */
1140  const char *type1 = cpl_propertylist_get_string(aPixtable->header, "CTYPE1"),
1141  *type2 = cpl_propertylist_get_string(aPixtable->header, "CTYPE2");
1142  cpl_ensure_code(type1 && type2 && !strncmp(type1, "RA---TAN", 9) &&
1143  !strncmp(type2, "DEC--TAN", 9),
1144  CPL_ERROR_UNSUPPORTED_MODE);
1145 
1146  cpl_msg_info(__func__, "Adapting WCS to RA/DEC=%f/%f deg", aRA, aDEC);
1147 
1148  /* delete the units of the x/y columns while we are working on them */
1149  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_XPOS, "");
1150  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_YPOS, "");
1151 
1152  /* replace coordinate values in the pixel table by the computed ones, *
1153  * carry out the spherical coordinate rotation for each pixel */
1154  float *xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
1155  *ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS);
1156  double dp = aDEC / CPL_MATH_DEG_RAD; /* delta_p in Paper II (in radians) */
1157  cpl_size n;
1158  #pragma omp parallel for default(none) /* as req. by Ralf */ \
1159  shared(aDEC, dp, nrow, xpos, ypos)
1160  for (n = 0; n < nrow; n++) {
1161  /* conversion from native spherical to celestial spherical coordinates *
1162  * alpha,delta as in Calabretta & Greisen 2002 A&A 395, 1077 (Paper II), *
1163  * Fig. 1; the formulae used are Eq. (2) in Sect. 2.3 in that paper and *
1164  * are generic but use that phi_p = 180 deg for zenithal projections *
1165  * (like TAN), i.e. cos(phi - phi_p) = -cos(phi) and similar for sin. */
1166  double phi = xpos[n],
1167  theta = ypos[n] + CPL_MATH_PI_2, /* add pi/2 again */
1168  ra = atan2(cos(theta) * sin(phi),
1169  sin(theta) * cos(dp) + cos(theta) * sin(dp) * cos(phi))
1170  * CPL_MATH_DEG_RAD,
1171  dec = asin(sin(theta) * sin(dp) - cos(theta) * cos(dp) * cos(phi))
1172  * CPL_MATH_DEG_RAD;
1173  /* The following should be *
1174  * xpos = aRA + ra; *
1175  * ypos = dec; *
1176  * but let's remove the zeropoint, to help accuracy despite *
1177  * the limited float precision (~7 digits, i.e. ~0.36 arcsec). */
1178  xpos[n] = ra;
1179  ypos[n] = dec - aDEC;
1180  } /* for n (table/matrix rows) */
1181 
1182  /* We produce output units in degrees (like wcslib); signify *
1183  * full WCS transformation by setting the "deg" output unit */
1184  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_XPOS, "deg");
1185  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_YPOS, "deg");
1186  /* append CRVAL now that we have a full WCS */
1187  cpl_propertylist_update_double(aPixtable->header, "CRVAL1", aRA);
1188  cpl_propertylist_update_double(aPixtable->header, "CRVAL2", aDEC);
1189  /* do the limit recomputation after setting the CRVALs */
1190  muse_pixtable_compute_limits(aPixtable);
1191 
1192  /* add the status header */
1193  cpl_propertylist_update_string(aPixtable->header, MUSE_HDR_PT_WCS,
1195  cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_WCS,
1196  MUSE_HDR_PT_WCS_COMMENT_POSI);
1197  return CPL_ERROR_NONE;
1198 } /* muse_wcs_position_celestial() */
1199 
1200 /*----------------------------------------------------------------------------*/
1223 /*----------------------------------------------------------------------------*/
1224 cpl_error_code
1225 muse_wcs_celestial_from_pixel(cpl_propertylist *aHeader, double aX, double aY,
1226  double *aRA, double *aDEC)
1227 {
1228  cpl_ensure_code(aHeader && aRA && aDEC, CPL_ERROR_NULL_INPUT);
1229  const char *type1 = cpl_propertylist_get_string(aHeader, "CTYPE1"),
1230  *type2 = cpl_propertylist_get_string(aHeader, "CTYPE2");
1231  cpl_ensure_code(type1 && type2 && !strncmp(type1, "RA---TAN", 9) &&
1232  !strncmp(type2, "DEC--TAN", 9),
1233  CPL_ERROR_UNSUPPORTED_MODE);
1234 
1235  muse_wcs *wcs = muse_wcs_new(aHeader);
1236  muse_wcs_celestial_from_pixel_fast(wcs, aX, aY, aRA, aDEC);
1237  cpl_free(wcs);
1238 
1239  return CPL_ERROR_NONE;
1240 } /* muse_wcs_celestial_from_pixel() */
1241 
1242 /*----------------------------------------------------------------------------*/
1267 /*----------------------------------------------------------------------------*/
1268 cpl_error_code
1269 muse_wcs_pixel_from_celestial(cpl_propertylist *aHeader, double aRA, double aDEC,
1270  double *aX, double *aY)
1271 {
1272  cpl_ensure_code(aHeader && aX && aY, CPL_ERROR_NULL_INPUT);
1273  /* make sure that the header represents TAN */
1274  const char *type1 = cpl_propertylist_get_string(aHeader, "CTYPE1"),
1275  *type2 = cpl_propertylist_get_string(aHeader, "CTYPE2");
1276  cpl_ensure_code(type1 && type2 && !strncmp(type1, "RA---TAN", 9) &&
1277  !strncmp(type2, "DEC--TAN", 9),
1278  CPL_ERROR_UNSUPPORTED_MODE);
1279 
1280  muse_wcs *wcs = muse_wcs_new(aHeader);
1281  if (wcs->cddet == 0.) { /* that's important here */
1282  *aX = *aY = NAN;
1283  cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
1284  cpl_free(wcs);
1285  return CPL_ERROR_ILLEGAL_INPUT;
1286  }
1287  wcs->crval1 /= CPL_MATH_DEG_RAD; /* convert to radians */
1288  wcs->crval2 /= CPL_MATH_DEG_RAD;
1289  muse_wcs_pixel_from_celestial_fast(wcs, aRA / CPL_MATH_DEG_RAD,
1290  aDEC / CPL_MATH_DEG_RAD, aX, aY);
1291  cpl_free(wcs);
1292 
1293  return CPL_ERROR_NONE;
1294 } /* muse_wcs_pixel_from_celestial() */
1295 
1296 /*----------------------------------------------------------------------------*/
1322 /*----------------------------------------------------------------------------*/
1323 cpl_error_code
1324 muse_wcs_projplane_from_celestial(cpl_propertylist *aHeader, double aRA,
1325  double aDEC, double *aX, double *aY)
1326 {
1327  cpl_ensure_code(aHeader && aX && aY, CPL_ERROR_NULL_INPUT);
1328  /* make sure that the header represents TAN */
1329  const char *type1 = cpl_propertylist_get_string(aHeader, "CTYPE1"),
1330  *type2 = cpl_propertylist_get_string(aHeader, "CTYPE2");
1331  cpl_ensure_code(type1 && type2 && !strncmp(type1, "RA---TAN", 9) &&
1332  !strncmp(type2, "DEC--TAN", 9),
1333  CPL_ERROR_UNSUPPORTED_MODE);
1334 
1335  /* spherical coordinate shift/translation */
1336  double a = aRA / CPL_MATH_DEG_RAD, /* RA in radians */
1337  d = aDEC / CPL_MATH_DEG_RAD, /* DEC in radians */
1338  /* alpha_p and delta_p in Paper II (in radians) */
1339  ap = cpl_propertylist_get_double(aHeader, "CRVAL1") / CPL_MATH_DEG_RAD,
1340  dp = cpl_propertylist_get_double(aHeader, "CRVAL2") / CPL_MATH_DEG_RAD,
1341  phi = atan2(-cos(d) * sin(a - ap),
1342  sin(d) * cos(dp) - cos(d) * sin(dp) * cos(a-ap))
1343  + 180 / CPL_MATH_DEG_RAD,
1344  theta = asin(sin(d) * sin(dp) + cos(d) * cos(dp) * cos(a-ap)),
1345  R_theta = CPL_MATH_DEG_RAD / tan(theta);
1346  /* spherical deprojection */
1347  *aX = R_theta * sin(phi),
1348  *aY = -R_theta * cos(phi);
1349 
1350  return CPL_ERROR_NONE;
1351 } /* muse_wcs_projplane_from_celestial() */
1352 
1353 /*----------------------------------------------------------------------------*/
1370 /*----------------------------------------------------------------------------*/
1371 cpl_error_code
1372 muse_wcs_projplane_from_pixel(cpl_propertylist *aHeader, double aX, double aY,
1373  double *aXOut, double *aYOut)
1374 {
1375  cpl_ensure_code(aHeader && aXOut && aYOut, CPL_ERROR_NULL_INPUT);
1376 
1377  muse_wcs *wcs = muse_wcs_new(aHeader);
1378  muse_wcs_projplane_from_pixel_fast(wcs, aX, aY, aXOut, aYOut);
1379  cpl_free(wcs);
1380 
1381  return CPL_ERROR_NONE;
1382 } /* muse_wcs_projplane_from_pixel() */
1383 
1384 /*----------------------------------------------------------------------------*/
1403 /*----------------------------------------------------------------------------*/
1404 cpl_error_code
1405 muse_wcs_pixel_from_projplane(cpl_propertylist *aHeader, double aX, double aY,
1406  double *aXOut, double *aYOut)
1407 {
1408  cpl_ensure_code(aHeader && aXOut && aYOut, CPL_ERROR_NULL_INPUT);
1409 
1410  muse_wcs *wcs = muse_wcs_new(aHeader);
1411  if (wcs->cddet == 0.) { /* that's important here */
1412  *aXOut = *aYOut = NAN;
1413  cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
1414  cpl_free(wcs);
1415  return CPL_ERROR_ILLEGAL_INPUT;
1416  }
1417  muse_wcs_pixel_from_projplane_fast(wcs, aX, aY, aXOut, aYOut);
1418  cpl_free(wcs);
1419 
1420  return CPL_ERROR_NONE;
1421 } /* muse_wcs_pixel_from_projplane() */
1422 
1423 /*----------------------------------------------------------------------------*/
1442 /*----------------------------------------------------------------------------*/
1443 cpl_error_code
1444 muse_wcs_get_angles(cpl_propertylist *aHeader, double *aXAngle, double *aYAngle)
1445 {
1446  cpl_ensure_code(aHeader && aXAngle && aYAngle, CPL_ERROR_NULL_INPUT);
1447 
1448  cpl_errorstate prestate = cpl_errorstate_get();
1449  double cd11 = cpl_propertylist_get_double(aHeader, "CD1_1"),
1450  cd22 = cpl_propertylist_get_double(aHeader, "CD2_2"),
1451  cd12 = cpl_propertylist_get_double(aHeader, "CD1_2"),
1452  cd21 = cpl_propertylist_get_double(aHeader, "CD2_1"),
1453  det = cd11 * cd22 - cd12 * cd21;
1454  cpl_ensure_code(cpl_errorstate_is_equal(prestate), cpl_error_get_code());
1455  if (det < 0.) {
1456  cd12 *= -1;
1457  cd11 *= -1;
1458  }
1459  if (cd12 == 0. && cd21 == 0.) { /* matrix without rotation */
1460  *aXAngle = 0.;
1461  *aYAngle = 0.;
1462  return CPL_ERROR_NONE;
1463  }
1464  /* angles in degrees */
1465  *aXAngle = atan2(cd12, cd11) * CPL_MATH_DEG_RAD;
1466  *aYAngle = atan2(-cd21, cd22) * CPL_MATH_DEG_RAD;
1467  return CPL_ERROR_NONE;
1468 } /* muse_wcs_get_angles() */
1469 
1470 /*----------------------------------------------------------------------------*/
1489 /*----------------------------------------------------------------------------*/
1490 cpl_error_code
1491 muse_wcs_get_scales(cpl_propertylist *aHeader, double *aXScale, double *aYScale)
1492 {
1493  cpl_ensure_code(aHeader && aXScale && aYScale, CPL_ERROR_NULL_INPUT);
1494 
1495  cpl_errorstate prestate = cpl_errorstate_get();
1496  double cd11 = cpl_propertylist_get_double(aHeader, "CD1_1"),
1497  cd22 = cpl_propertylist_get_double(aHeader, "CD2_2"),
1498  cd12 = cpl_propertylist_get_double(aHeader, "CD1_2"),
1499  cd21 = cpl_propertylist_get_double(aHeader, "CD2_1"),
1500  det = cd11 * cd22 - cd12 * cd21;
1501  cpl_ensure_code(cpl_errorstate_is_equal(prestate), cpl_error_get_code());
1502 
1503  if (det < 0.) {
1504  cd12 *= -1;
1505  cd11 *= -1;
1506  }
1507  if (cd12 == 0. && cd21 == 0.) { /* matrix without rotation */
1508  *aXScale = cd11;
1509  *aYScale = cd22;
1510  return CPL_ERROR_NONE;
1511  }
1512  *aXScale = sqrt(cd11*cd11 + cd12*cd12); /* only the absolute value */
1513  *aYScale = sqrt(cd22*cd22 + cd21*cd21);
1514  return CPL_ERROR_NONE;
1515 } /* muse_wcs_get_scales() */
1516 
1517 /*----------------------------------------------------------------------------*/
1534 /*----------------------------------------------------------------------------*/
1535 muse_wcs *
1536 muse_wcs_new(cpl_propertylist *aHeader)
1537 {
1538  muse_wcs *wcs = cpl_calloc(1, sizeof(muse_wcs));
1539  if (!aHeader) {
1540  wcs->cd11 = wcs->cd22 = wcs->cddet = 1.; /* see below */
1541  return wcs;
1542  }
1543 
1544  cpl_errorstate prestate = cpl_errorstate_get();
1545  wcs->crpix1 = cpl_propertylist_get_double(aHeader, "CRPIX1");
1546  wcs->crpix2 = cpl_propertylist_get_double(aHeader, "CRPIX2");
1547  wcs->crval1 = cpl_propertylist_get_double(aHeader, "CRVAL1");
1548  wcs->crval2 = cpl_propertylist_get_double(aHeader, "CRVAL2");
1549  if (!cpl_errorstate_is_equal(prestate)) {
1550  /* all these headers default to 0.0 following the FITS *
1551  * Standard, so we can ignore any errors set here */
1552  cpl_errorstate_set(prestate);
1553  }
1554 
1555  prestate = cpl_errorstate_get();
1556  wcs->cd11 = cpl_propertylist_get_double(aHeader, "CD1_1");
1557  wcs->cd22 = cpl_propertylist_get_double(aHeader, "CD2_2");
1558  wcs->cd12 = cpl_propertylist_get_double(aHeader, "CD1_2");
1559  wcs->cd21 = cpl_propertylist_get_double(aHeader, "CD2_1");
1560  if (!cpl_errorstate_is_equal(prestate) &&
1561  wcs->cd11 == 0. && wcs->cd12 == 0. && wcs->cd21 == 0. && wcs->cd22 == 0.) {
1562  /* FITS Standard says to handle the CD matrix like the PC *
1563  * matrix in this case, with 1 for the diagonal elements */
1564  wcs->cd11 = wcs->cd22 = wcs->cddet = 1.;
1565  cpl_errorstate_set(prestate); /* not a real error */
1566  }
1567  wcs->cddet = wcs->cd11 * wcs->cd22 - wcs->cd12 * wcs->cd21;
1568  if (wcs->cddet == 0.) {
1569  cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
1570  }
1571 
1572  /* wcs->iscelsph defaults to 0 = CPL_FALSE, leave it at that */
1573  if (getenv("MUSE_DEBUG_WCS") && atoi(getenv("MUSE_DEBUG_WCS")) > 0) {
1574  cpl_msg_debug(__func__, "wcs: axis1 = { %f %f %e }, axis2 = { %f %f %e }, "
1575  "crossterms = { %e %e }, det = %e",
1576  wcs->crpix1, wcs->crval1, wcs->cd11,
1577  wcs->crpix2, wcs->crval2, wcs->cd22,
1578  wcs->cd12, wcs->cd21, wcs->cddet);
1579  }
1580  return wcs;
1581 } /* muse_wcs_new() */
1582 
double crpix2
Definition: muse_wcs.h:94
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:1491
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:1405
void muse_image_delete(muse_image *aImage)
Deallocate memory associated to a muse_image object.
Definition: muse_image.c:84
#define MUSE_HDR_PT_XLO
FITS header keyword contains the lower limit of the data in x-direction.
cpl_size muse_pixtable_get_nrow(muse_pixtable *aPixtable)
get the number of rows within the pixel table
A structure containing a spatial two-axis WCS.
Definition: muse_wcs.h:93
double muse_pfits_get_ra(const cpl_propertylist *aHeaders)
find out the right ascension
Definition: muse_pfits.c:206
double cd22
Definition: muse_wcs.h:96
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:47
muse_image * muse_combine_median_create(muse_imagelist *aImages)
Median combine a list of input images.
Definition: muse_combine.c:317
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:233
cpl_image * stat
the statistics extension
Definition: muse_image.h:65
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:1324
#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:58
Structure definition of MUSE three extension FITS file.
Definition: muse_image.h:41
cpl_table * table
The pixel table.
cpl_propertylist * header
the FITS header
Definition: muse_image.h:73
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:992
muse_resampling_crstats_type crtype
cpl_image * dq
the data quality extension
Definition: muse_image.h:57
cpl_table * muse_cpltable_new(const muse_cpltable_def *aDef, cpl_size aLength)
Create an empty table according to the specified definition.
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:1536
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:858
muse_wcs_centroid_type
Type of centroiding algorithm to use.
Definition: muse_wcs.h:78
#define MUSE_HDR_PT_PREDAR_YHI
Structure definition of MUSE pixel table.
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:1225
#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:1473
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:164
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:126
double cddet
Definition: muse_wcs.h:97
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:582
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:1372
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:840
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:1269
double muse_pfits_get_dec(const cpl_propertylist *aHeaders)
find out the declination
Definition: muse_pfits.c:224
#define MUSE_HDR_PT_PREDAR_XHI
muse_datacube * muse_resampling_cube(muse_pixtable *aPixtable, muse_resampling_params *aParams, muse_pixgrid **aPixgrid)
Resample a pixel table onto a regular grid structure representing a FITS NAXIS=3 datacube.
cpl_imagelist * dq
the optional cube containing the bad pixel status
Definition: muse_datacube.h:81
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:404
#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:396
cpl_error_code muse_image_reject_from_dq(muse_image *aImage)
Reject pixels of a muse_image depending on its DQ data.
Definition: muse_image.c:848
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:204
#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:903
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:829
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:1444
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:65
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:1132
double crval2
Definition: muse_wcs.h:95
#define MUSE_HDR_PT_WCS
#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.
muse_ins_mode muse_pfits_get_mode(const cpl_propertylist *aHeaders)
find out the observation mode
Definition: muse_pfits.c:1003
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:198
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:838