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