MUSE Pipeline Reference Manual  2.1.1
muse_postproc.c
1 /* -*- Mode: C; tab-width: 2; indent-tabs-mode: nil; c-basic-offset: 2 -*- */
2 /* vim:set sw=2 sts=2 et cin: */
3 /*
4  * This file is part of the MUSE Instrument Pipeline
5  * Copyright (C) 2005-2015 European Southern Observatory
6  *
7  * This program is free software; you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation; either version 2 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this program; if not, write to the Free Software
19  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
20  * 02110-1301, USA.
21  */
22 
23 #ifdef HAVE_CONFIG_H
24 #include <config.h>
25 #endif
26 
27 /*----------------------------------------------------------------------------*
28  * Includes *
29  *----------------------------------------------------------------------------*/
30 #include <string.h>
31 #include <strings.h>
32 #include "muse_postproc.h"
33 #include "muse_instrument.h"
34 
35 #include "muse_astro.h"
36 #include "muse_dar.h"
37 #include "muse_flux.h"
38 #include "muse_imagelist.h"
39 #include "muse_pfits.h"
40 #include "muse_phys.h"
41 #include "muse_rvcorrect.h"
42 #include "muse_utils.h"
43 #include "muse_wcs.h"
44 #include "muse_idp.h"
45 
46 #include "muse_lsf_params.h"
47 #ifdef USE_LSF_PARAMS
48 #include "muse_sky_old.h"
49 #endif
50 
51 #undef CREATE_IDP
52 
53 /*----------------------------------------------------------------------------*/
60 /*----------------------------------------------------------------------------*/
61 
64 /*----------------------------------------------------------------------------*/
75 /*----------------------------------------------------------------------------*/
78 {
79  muse_postproc_properties *prop = cpl_calloc(1, sizeof(muse_postproc_properties));
80  switch (aType) {
81  case MUSE_POSTPROC_SCIPOST:
82  case MUSE_POSTPROC_STANDARD:
83  case MUSE_POSTPROC_ASTROMETRY:
84  prop->type = aType;
85  break;
86  default:
87  cpl_msg_error(__func__, "No such setup known: %d", aType);
88  cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
89  cpl_free(prop);
90  return NULL;
91  }
92  /* rvtype of 0 is already MUSE_RVCORRECT_NONE */
93  return prop;
94 } /* muse_postproc_properties_new() */
95 
96 /*----------------------------------------------------------------------------*/
102 /*----------------------------------------------------------------------------*/
103 void
105 {
106  if (!aProp) {
107  return;
108  }
109  cpl_table_delete(aProp->exposures);
110  muse_table_delete(aProp->response);
111  muse_table_delete(aProp->telluric);
112  cpl_table_delete(aProp->extinction);
113  cpl_propertylist_delete(aProp->wcs);
115 #ifdef USE_LSF_PARAMS
117 #endif
118  cpl_table_delete(aProp->sky_lines);
119  cpl_table_delete(aProp->sky_continuum);
120  muse_mask_delete(aProp->sky_mask);
121  cpl_frame_delete(aProp->refframe);
122  cpl_table_delete(aProp->tellregions);
123  cpl_free(aProp);
124 } /* muse_postproc_properties_delete() */
125 
126 /*----------------------------------------------------------------------------*/
133 /*----------------------------------------------------------------------------*/
134 cpl_boolean
135 muse_postproc_check_save_param(const char *aSave, const char *aValid)
136 {
137  cpl_ensure(aSave, CPL_ERROR_NULL_INPUT, CPL_FALSE);
138  if (strlen(aSave) < 4) { /* less than "cube"! */
139  cpl_error_set_message(__func__, CPL_ERROR_DATA_NOT_FOUND,
140  "no (valid) save option given!");
141  return CPL_FALSE;
142  }
143  /* now compare all given output products against the allowed options */
144  cpl_boolean allvalid = CPL_TRUE;
145  cpl_array *asave = muse_cplarray_new_from_delimited_string(aSave, ","),
146  *avalid = muse_cplarray_new_from_delimited_string(aValid, ",");
147  int i, j,
148  ns = cpl_array_get_size(asave),
149  nv = cpl_array_get_size(avalid);
150  for (i = 0; i < ns; i++) {
151  cpl_boolean valid = CPL_FALSE;
152  for (j = 0; j < nv; j++) { /* go through all valid options */
153  if (!strcmp(cpl_array_get_string(asave, i),
154  cpl_array_get_string(avalid, j))) {
155  valid = CPL_TRUE;
156  }
157  } /* for j (all valid strings) */
158  if (!valid) { /* this string was apparently not valid! */
159  cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
160  "save option %d (%s) is not valid!", i + 1,
161  cpl_array_get_string(asave, i));
162  allvalid = CPL_FALSE;
163  }
164  } /* for i (all given strings) */
165  cpl_array_delete(asave);
166  cpl_array_delete(avalid);
167  return allvalid;
168 } /* muse_scipost_check_save_param() */
169 
170 /*----------------------------------------------------------------------------*/
185 /*----------------------------------------------------------------------------*/
187 muse_postproc_get_resampling_type(const char *aResampleString)
188 {
189  cpl_ensure(aResampleString, CPL_ERROR_NULL_INPUT, MUSE_RESAMPLE_WEIGHTED_DRIZZLE);
190  if (!strncmp(aResampleString, "nearest", 8)) {
191  return MUSE_RESAMPLE_NEAREST;
192  }
193  if (!strncmp(aResampleString, "linear", 7)) {
195  }
196  if (!strncmp(aResampleString, "quadratic", 10)) {
198  }
199  if (!strncmp(aResampleString, "renka", 6)) {
201  }
202  if (!strncmp(aResampleString, "drizzle", 8)) {
204  }
205  if (!strncmp(aResampleString, "lanczos", 8)) {
207  }
208  return MUSE_RESAMPLE_WEIGHTED_DRIZZLE; /* catch all fallback */
209 } /* muse_postproc_get_resampling_type() */
210 
211 /*----------------------------------------------------------------------------*/
226 /*----------------------------------------------------------------------------*/
228 muse_postproc_get_cr_type(const char *aCRTypeString)
229 {
230  cpl_ensure(aCRTypeString, CPL_ERROR_NULL_INPUT, MUSE_RESAMPLING_CRSTATS_IRAF);
231  if (!strncmp(aCRTypeString, "iraf", 5)) {
233  }
234  if (!strncmp(aCRTypeString, "mean", 5)) {
236  }
237  if (!strncmp(aCRTypeString, "median", 7)) {
239  }
240  return MUSE_RESAMPLING_CRSTATS_MEDIAN; /* catch all fallback */
241 } /* muse_postproc_get_cr_type() */
242 
243 /*----------------------------------------------------------------------------*/
257 /*----------------------------------------------------------------------------*/
259 muse_postproc_get_cube_format(const char *aFormatString)
260 {
261  cpl_ensure(aFormatString, CPL_ERROR_NULL_INPUT, MUSE_CUBE_TYPE_FITS);
262  if (!strncmp(aFormatString, "Cube", 5)) {
263  return MUSE_CUBE_TYPE_FITS;
264  }
265  if (!strncmp(aFormatString, "Euro3D", 7)) {
266  return MUSE_CUBE_TYPE_EURO3D;
267  }
268  if (!strncmp(aFormatString, "xCube", 6)) {
269  return MUSE_CUBE_TYPE_FITS_X;
270  }
271  if (!strncmp(aFormatString, "xEuro3D", 8)) {
273  }
274  if (!strncmp(aFormatString, "sdpCube", 8)) {
275  return MUSE_CUBE_TYPE_SDP;
276  }
277  return MUSE_CUBE_TYPE_FITS; /* catch all fallback */
278 } /* muse_postproc_get_cube_format() */
279 
280 /*----------------------------------------------------------------------------*/
294 /*----------------------------------------------------------------------------*/
296 muse_postproc_get_weight_type(const char *aWeightString)
297 {
298  cpl_ensure(aWeightString, CPL_ERROR_NULL_INPUT, MUSE_XCOMBINE_EXPTIME);
299  if (!strncmp(aWeightString, "exptime", 8)) {
300  return MUSE_XCOMBINE_EXPTIME;
301  }
302  if (!strncmp(aWeightString, "fwhm", 5)) {
303  return MUSE_XCOMBINE_FWHM;
304  }
305  if (!strncmp(aWeightString, "none", 5)) {
306  return MUSE_XCOMBINE_NONE;
307  }
308  return MUSE_XCOMBINE_EXPTIME; /* catch all fallback */
309 } /* muse_postproc_get_weight_type() */
310 
311 /*----------------------------------------------------------------------------*/
336 /*----------------------------------------------------------------------------*/
337 cpl_table *
338 muse_postproc_load_nearest(const cpl_propertylist *aHeader,
339  const cpl_frame *aFrame, float aWarnLimit,
340  float aErrLimit, double *aRA, double *aDEC)
341 {
342  cpl_ensure(aHeader && aFrame, CPL_ERROR_NULL_INPUT, NULL);
343  cpl_errorstate state = cpl_errorstate_get();
344  double raref = muse_pfits_get_ra(aHeader),
345  decref = muse_pfits_get_dec(aHeader);
346  cpl_ensure(cpl_errorstate_is_equal(state), CPL_ERROR_DATA_NOT_FOUND, NULL);
347  cpl_msg_debug(__func__, "reference coordinates: RA = %e deg, DEC =%e deg",
348  raref, decref);
349  if (aRA) {
350  *aRA = raref;
351  }
352  if (aDEC) {
353  *aDEC = decref;
354  }
355 
356  const char *fn = cpl_frame_get_filename(aFrame);
357  cpl_propertylist *header;
358  double dmin = FLT_MAX; /* minimum distance yet found in deg */
359  int iext, inearest = -1, next = cpl_fits_count_extensions(fn);
360  for (iext = 1; iext <= next; iext++) {
361  header = cpl_propertylist_load(fn, iext);
362  const char *extname = cpl_propertylist_get_string(header, "EXTNAME");
363  double ra = muse_pfits_get_ra(header),
364  dec = muse_pfits_get_dec(header),
365  d = muse_astro_angular_distance(ra, dec, raref, decref);
366  cpl_msg_debug(__func__, "extension %d [%s]: RA = %e deg, DEC = %e deg --> "
367  "d = %e deg", iext, extname, ra, dec, d);
368  if (d < dmin) {
369  dmin = d;
370  inearest = iext;
371  }
372  cpl_propertylist_delete(header);
373  } /* for iext */
374  /* warn for distances larger than x arcsec */
375  if (dmin * 3600. > aErrLimit) {
376  char *msg = cpl_sprintf("Distance of nearest reference table to observed "
377  "position is %.2f arcmin, certainly a wrong "
378  "reference object!", dmin * 60.);
379  cpl_msg_error(__func__, "%s", msg);
380  cpl_error_set_message(__func__, CPL_ERROR_INCOMPATIBLE_INPUT, "%s", msg);
381  cpl_free(msg);
382  return NULL;
383  }
384  if (dmin * 3600. > aWarnLimit) {
385  cpl_msg_warning(__func__, "Distance of nearest reference table to observed "
386  "position is %.2f arcmin! (Wrong reference object?)",
387  dmin * 60.);
388  }
389  header = cpl_propertylist_load(fn, inearest);
390  const char *extname = cpl_propertylist_get_string(header, "EXTNAME");
391  cpl_msg_info(__func__, "Loading \"%s[%s]\" (distance %.1f arcsec)", fn,
392  extname, dmin * 3600.);
393  cpl_propertylist_delete(header);
394  cpl_table *table = cpl_table_load(fn, inearest, 1);
395  return table;
396 } /* muse_postproc_load_nearest() */
397 
398 /*----------------------------------------------------------------------------*/
414 /*----------------------------------------------------------------------------*/
415 cpl_error_code
417  const muse_table *aResponse)
418 {
419  cpl_ensure_code(aPt && aPt->header, CPL_ERROR_NULL_INPUT);
420  if (!aResponse) {
421  return CPL_ERROR_NONE;
422  }
423  /* cross-check the header of the response curve against *
424  * the flat-field spectrum status of the pixel table */
425  cpl_boolean respff = cpl_propertylist_has(aResponse->header,
426  MUSE_HDR_FLUX_FFCORR),
427  ptff = cpl_propertylist_has(aPt->header, MUSE_HDR_PT_FFCORR);
428  if (respff == ptff || !ptff) {
429  /* same state for both means we don't need to do anything, same if *
430  * the pixel table was not corrected for ffspec */
431  return CPL_ERROR_NONE;
432  }
433 
434  cpl_msg_warning(__func__, "Adapt pixel table to %s for backward compatibility"
435  ": revert correction by flat-field spectrum!",
436  MUSE_TAG_STD_RESPONSE);
437 
438  /* so we need to divide again by the flat-field spectrum */
439  cpl_array *lambdas = muse_cpltable_extract_column(aPt->ffspec,
441  *ffdata = muse_cpltable_extract_column(aPt->ffspec,
443  muse_pixtable_spectrum_apply(aPt, lambdas, ffdata,
444  MUSE_PIXTABLE_OPERATION_DIVIDE);
445  cpl_array_unwrap(lambdas);
446  cpl_array_unwrap(ffdata);
447  cpl_table_delete(aPt->ffspec);
448  aPt->ffspec = NULL;
449  cpl_propertylist_erase(aPt->header, MUSE_HDR_PT_FFCORR);
450  cpl_msg_info(__func__, "Pixel table now convolved with flat-field spectrum "
451  "again, removed %s keyword from header.", MUSE_HDR_PT_FFCORR);
452  return CPL_ERROR_NONE;
453 } /* muse_postproc_revert_ffspec_maybe() */
454 
455 /*----------------------------------------------------------------------------*/
498 /*----------------------------------------------------------------------------*/
499 void *
501  unsigned int aIndex,
502  muse_postproc_sky_outputs *aSkyOut)
503 {
504  cpl_ensure(aProp && aProp->exposures, CPL_ERROR_NULL_INPUT, NULL);
505  cpl_ensure((int)aIndex < cpl_table_get_nrow(aProp->exposures),
506  CPL_ERROR_ILLEGAL_INPUT, NULL);
507  cpl_msg_info(__func__, "Starting to process exposure %d (DATE-OBS=%s)",
508  aIndex + 1, cpl_table_get_string(aProp->exposures, "DATE-OBS",
509  aIndex));
510  cpl_ensure(aProp->exposures, CPL_ERROR_NULL_INPUT, NULL);
511 
512  cpl_table *exposure = cpl_table_extract(aProp->exposures, aIndex, 1);
514  aProp->lambdamin,
515  aProp->lambdamax);
516  cpl_table_delete(exposure);
517  if (!pt) {
518  return NULL;
519  }
520  /* erase pre-existing QC parameters */
521  cpl_propertylist_erase_regexp(pt->header, "ESO QC ", 0);
522 
523  /* see if we need to revert the flat-field spectrum correction, *
524  * to match pixel table and response curve */
526 
527  /* try to recognize a lab exposure by missing RA/DEC or weird airmass */
528  cpl_boolean labdata = CPL_FALSE;
529  cpl_errorstate prestate = cpl_errorstate_get();
530  double airmass = muse_astro_airmass(pt->header);
531  double ra = muse_pfits_get_ra(pt->header),
532  dec = muse_pfits_get_dec(pt->header);
533  if (!cpl_errorstate_is_equal(prestate) || airmass < 1.) {
534  cpl_msg_warning(__func__, "This seems to be lab data (RA=%f DEC=%f, "
535  "airmass=%f), not doing any DAR correction!", ra, dec,
536  airmass);
537  cpl_errorstate_dump(prestate, CPL_FALSE, NULL);
538  cpl_errorstate_set(prestate);
539  labdata = CPL_TRUE;
540  }
541  cpl_error_code rc = CPL_ERROR_NONE;
542  if (!labdata && muse_pfits_get_mode(pt->header) <= MUSE_MODE_WFM_AO_N) {
543  cpl_msg_debug(__func__, "WFM detected: starting DAR correction");
544  rc = muse_dar_correct(pt, aProp->lambdaref);
545  cpl_msg_debug(__func__, "DAR correction returned rc=%d: %s", rc,
546  rc != CPL_ERROR_NONE ? cpl_error_get_message() : "");
547  /* check and possibly correct the DAR quality, if requested to do so */
548  if (aProp->darcheck != MUSE_POSTPROC_DARCHECK_NONE) {
549  cpl_boolean dorefine = aProp->darcheck == MUSE_POSTPROC_DARCHECK_CORRECT;
550  double maxshift = 0;
551  rc = muse_dar_check(pt, &maxshift, dorefine, NULL);
552  if (rc != CPL_ERROR_NONE) {
553  cpl_msg_warning(__func__, "Maximum %s shift %f\'\' (rc = %d: %s)",
554  dorefine ? "corrected" : "detected", maxshift, rc,
555  cpl_error_get_message());
556  } else {
557  cpl_msg_info(__func__, "Maximum %s shift %f\'\'",
558  dorefine ? "corrected" : "detected", maxshift);
559  }
560  } /* if DAR checking */
561  } /* if not labdata but WFM */
562 
563  /* do flux calibration (before sky subtraction!), if possible */
564  if (aProp->type != MUSE_POSTPROC_STANDARD) {
565  if (aProp->response) {
566  rc = muse_flux_calibrate(pt, aProp->response, aProp->extinction,
567  aProp->telluric);
568  if (rc != CPL_ERROR_NONE && rc != CPL_ERROR_UNSUPPORTED_MODE) {
569  /* for a more serious error, output a warning at least *
570  * (the data may still be usable) */
571  cpl_msg_warning(__func__, "Flux calibration returned rc=%d: %s", rc,
572  cpl_error_get_message());
573  } else {
574  /* for no error or unsupported mode (an already flux *
575  * calibrated table!) a debug message seems more appropriate */
576  cpl_msg_debug(__func__, "Flux calibration returned rc=%d: %s", rc,
577  rc != CPL_ERROR_NONE ? cpl_error_get_message_default(rc) : "");
578  }
579  } else {
580  cpl_msg_info(__func__, "Skipping flux calibration (no %s curve)",
581  MUSE_TAG_STD_RESPONSE);
582  } /* else */
583  } /* if not POSTPROC_STANDARD */
584 
585  cpl_table *sky_cont = NULL,
586  *sky_lines = NULL;
587  if ((aProp->skymethod == MUSE_POSTPROC_SKYMETHOD_MODEL && aProp->sky_lines) ||
588  (aProp->skymethod == MUSE_POSTPROC_SKYMETHOD_SIMPLE)) {
589  /* Process optional sky mask */
591  if (aProp->sky_mask) {
592  cpl_table_select_all(sky_pt->table);
594  }
595  /* Create white light image and sky mask */
596  cpl_msg_info(__func__, "Intermediate resampling to produce white-light image"
597  " for sky-spectrum creation:");
598  cpl_msg_indent_more();
601  params->crsigma = aProp->skymodel_params.crsigmac;
602  muse_datacube *cube = muse_resampling_cube(sky_pt, params, NULL);
603  if (!cube) {
604  cpl_msg_error(__func__, "Could not create cube for whitelight image");
605  muse_pixtable_delete(sky_pt);
607  return NULL;
608  }
609  muse_table *fwhite = muse_table_load_filter(NULL, "white");
610  muse_image *whitelight = muse_datacube_collapse(cube, fwhite);
611  muse_table_delete(fwhite);
612  cpl_msg_indent_less();
613  muse_mask *sky_mask = muse_sky_create_skymask(whitelight,
614  aProp->skymodel_params.ignore,
615  aProp->skymodel_params.fraction,
616  /* hardcode only sensible prefix */ "ESO QC SCIPOST");
617 
618  /* Apply mask and create spectrum, possibly with one CR-reject iteration */
619  cpl_table_select_all(sky_pt->table);
620  muse_pixtable_and_selected_mask(sky_pt, sky_mask);
621  cpl_table *spectrum =
622  muse_resampling_spectrum_iterate(sky_pt, aProp->skymodel_params.sampling,
623  0., aProp->skymodel_params.crsigmas, 1);
624  if (aProp->skymethod == MUSE_POSTPROC_SKYMETHOD_MODEL) {
625  /* Fit sky lines spectrum */
626  cpl_array *lambda = muse_cpltable_extract_column(spectrum, "lambda");
627  double lambda_low = cpl_array_get_min(lambda);
628  double lambda_high = cpl_array_get_max(lambda);
629  sky_lines = cpl_table_duplicate(aProp->sky_lines);
630  muse_sky_lines_set_range(sky_lines, lambda_low-5, lambda_high+5);
631 
632  // do the fit, ignoring possible errors
633  prestate = cpl_errorstate_get();
634  if (aProp->lsf_cube) {
635  cpl_image *lsfImage = muse_lsf_average_cube_all(aProp->lsf_cube, sky_pt);
636  muse_wcs *lsfWCS = muse_lsf_cube_get_wcs_all(aProp->lsf_cube);
637  /* Filter the image with a rectangle to model the spectrum binning */
638  muse_lsf_fold_rectangle(lsfImage, lsfWCS, aProp->skymodel_params.sampling);
639 
640  muse_sky_lines_fit(spectrum, sky_lines, lsfImage, lsfWCS);
641  if (!cpl_errorstate_is_equal(prestate)) {
642  cpl_errorstate_dump(prestate, CPL_FALSE, NULL);
643  cpl_errorstate_set(prestate);
644  }
645 
646  /* Create continuum */
647  if (!aProp->sky_continuum) {
648  cpl_msg_info(__func__, "No sky continuum given, create it from the data");
649  sky_cont = muse_sky_continuum_create(spectrum, sky_lines, lsfImage, lsfWCS,
650  aProp->skymodel_params.csampling);
651  } else {
652  cpl_msg_info(__func__, "Using sky continuum given as input");
653  }
654 
655  /* clean up */
656  cpl_image_delete(lsfImage);
657 #ifdef USE_LSF_PARAMS
658  } else if (aProp->lsf_params) { // Old LSF params code
659  // do the fit, ignoring possible errors
660  muse_sky_lines_fit_old(spectrum, sky_lines);
661  if (!cpl_errorstate_is_equal(prestate)) {
662  cpl_errorstate_dump(prestate, CPL_FALSE, NULL);
663  cpl_errorstate_set(prestate);
664  }
665 
666  if (!aProp->sky_continuum) {
667  cpl_msg_info(__func__, "No sky continuum given, create it from the data");
668  muse_sky_subtract_lines_old(sky_pt, sky_lines, aProp->lsf_params);
669  sky_cont = muse_resampling_spectrum(sky_pt,
670  aProp->skymodel_params.csampling);
671  cpl_table_erase_column(sky_cont, "stat");
672  cpl_table_erase_column(sky_cont, "dq");
673  cpl_table_name_column(sky_cont, "data", "flux");
674  }
675 #endif
676  }
677  cpl_array_unwrap(lambda);
678  } else {
679  /* simple subtraction: pretend that the spectrum is a *
680  * continuum and use that function for direct subtraction */
681  cpl_table *spec2 = cpl_table_duplicate(spectrum);
682  cpl_table_name_column(spec2, "data", "flux");
683  rc = muse_sky_subtract_continuum(pt, spec2);
684  cpl_table_delete(spec2);
685  } /* else: simple */
686  if (aSkyOut) {
687  aSkyOut->mask = sky_mask;
688  aSkyOut->spectrum = spectrum;
689  } else {
690  muse_mask_delete(sky_mask);
691  cpl_table_delete(spectrum);
692  }
693  muse_image_delete(whitelight);
695  muse_datacube_delete(cube);
696  muse_pixtable_delete(sky_pt);
697  } /* if MUSE_POSTPROC_SKYMETHOD_MODEL or MUSE_POSTPROC_SKYMETHOD_SIMPLE */
698  /* In the following we do the sky subtraction, if we have the products *
699  * for it. The SKYMETHOD_NONE here reflects the subtract-model option *
700  * of the muse_scipost recipe. */
701  if ((aProp->skymethod == MUSE_POSTPROC_SKYMETHOD_NONE) ||
702  (aProp->skymethod == MUSE_POSTPROC_SKYMETHOD_MODEL)) {
703  if (aProp->sky_continuum) {
704  cpl_msg_debug(__func__, "doing sky continuum subtraction with input "
705  "spectrum");
707  } else if (sky_cont) {
708  cpl_msg_debug(__func__, "doing sky continuum subtraction with created "
709  "spectrum");
710  rc = muse_sky_subtract_continuum(pt, sky_cont);
711  }
712  if (sky_lines) {
713  cpl_msg_debug(__func__, "doing sky line subtraction with created list");
714  if (aProp->lsf_cube != NULL) {
715  rc = muse_sky_subtract_lines(pt, sky_lines, aProp->lsf_cube);
716 #ifdef USE_LSF_PARAMS
717  } else if (aProp->lsf_params != NULL) {
718  rc = muse_sky_subtract_lines_old(pt, sky_lines, aProp->lsf_params);
719 #endif
720  }
721  } else if (aProp->sky_lines) {
722  cpl_msg_debug(__func__, "doing sky line subtraction with input list");
723  if (aProp->lsf_cube != NULL) {
724  rc = muse_sky_subtract_lines(pt, aProp->sky_lines, aProp->lsf_cube);
725 #ifdef USE_LSF_PARAMS
726  } else if (aProp->lsf_params != NULL) {
727  rc = muse_sky_subtract_lines_old(pt, aProp->sky_lines, aProp->lsf_params);
728 #endif
729  }
730  }
731  } /* if sky modeling */
732  /* either transfer the sky lines and continuum out or clean *
733  * them up; both should still work, if one or both are NULL */
734  if (aSkyOut) {
735  aSkyOut->lines = sky_lines;
736  aSkyOut->continuum = sky_cont;
737  } else {
738  cpl_table_delete(sky_lines);
739  cpl_table_delete(sky_cont);
740  }
741 
742  /* carry out radial velocity correction, if specified */
743  muse_rvcorrect(pt, aProp->rvtype);
744 
745  /* do flux response computation and return */
746  if (aProp->type == MUSE_POSTPROC_STANDARD) {
747  double raref, decref;
748  cpl_table *reftable = muse_postproc_load_nearest(pt->header,
749  aProp->refframe, 20., 60.,
750  &raref, &decref);
751  if (!reftable) {
752  cpl_msg_error(__func__, "Correct %s could not be loaded, observed "
753  "target not listed?", MUSE_TAG_STD_FLUX_TABLE);
755  return NULL;
756  }
757  if (muse_flux_reference_table_check(reftable) != CPL_ERROR_NONE) {
758  cpl_msg_error(__func__, "%s does not have a recognized format!",
759  MUSE_TAG_STD_FLUX_TABLE);
760  cpl_table_delete(reftable);
762  return NULL;
763  }
764 
765  prestate = cpl_errorstate_get();
767  rc = muse_flux_integrate_std(pt, aProp->profile, flux);
768  if (flux->intimage) {
769  cpl_msg_debug(__func__, "Flux integration found %"CPL_SIZE_FORMAT" stars",
770  cpl_image_get_size_y(flux->intimage->data));
771  }
773  if (airmass < 0) {
774  cpl_msg_warning(__func__, "Invalid airmass derived (%.3e), resetting to "
775  " zero", airmass);
776  airmass = 0.; /* set to zero to leave out uncertain extinction term */
777  }
778  flux->raref = raref;
779  flux->decref = decref;
780  rc = muse_flux_response_compute(flux, aProp->select, airmass, reftable,
781  aProp->tellregions, aProp->extinction);
782  cpl_table_delete(reftable);
783  rc = muse_flux_get_response_table(flux, aProp->smooth);
784  rc = muse_flux_get_telluric_table(flux);
785  rc = muse_flux_compute_qc(flux);
786 
787  if (!cpl_errorstate_is_equal(prestate)) {
788  cpl_msg_warning(__func__, "The following errors happened while computing "
789  "the flux calibration (rc = %d/%s):", rc,
790  cpl_error_get_message_default(rc));
791  cpl_errorstate_dump(prestate, CPL_FALSE, muse_cplerrorstate_dump_some);
792  }
793  return flux;
794  } /* if MUSE_POSTPROC_STANDARD */
795 
796  /* do astrometry computation and return */
797  if (aProp->type == MUSE_POSTPROC_ASTROMETRY) {
798  cpl_table *reftable = muse_postproc_load_nearest(pt->header,
799  aProp->refframe, 20., 60.,
800  NULL, NULL);
801  if (!reftable) {
802  cpl_msg_error(__func__, "Correct %s could not be loaded, observed "
803  "target not listed?", MUSE_TAG_ASTROMETRY_REFERENCE);
805  return NULL;
806  }
807  prestate = cpl_errorstate_get();
808 
810  wcs->crpix1 = aProp->crpix1; /* set center of rotation */
811  wcs->crpix2 = aProp->crpix2;
812  rc = muse_wcs_locate_sources(pt, fabsf(aProp->detsigma), aProp->centroid,
813  wcs);
814  if (aProp->detsigma < 0.) {
815  /* loop through possible detection sigmas to find the optimal solution */
816  rc = muse_wcs_optimize_solution(wcs, aProp->detsigma, aProp->centroid,
817  reftable, aProp->radius, aProp->faccuracy,
818  aProp->niter, aProp->rejsigma);
819  } else {
820  /* use just the given detection sigma for the solution */
821  rc = muse_wcs_solve(wcs, reftable, aProp->radius, aProp->faccuracy,
822  aProp->niter, aProp->rejsigma);
823  if (wcs->wcs) { /* set the keyword that muse_wcs_solve() does not save */
824  cpl_propertylist_update_float(wcs->wcs, MUSE_HDR_WCS_DETSIGMA,
825  aProp->detsigma);
826  cpl_propertylist_set_comment(wcs->wcs, MUSE_HDR_WCS_DETSIGMA,
827  MUSE_HDR_WCS_DETSIGMA_C_ONE);
828  }
829  }
830  cpl_table_delete(reftable);
831  cpl_msg_debug(__func__, "Solving astrometry (initial radius %.1f pix, "
832  "initial relative accuracy %.1f, %d iterations with rejection"
833  " sigma %.2f) returned rc=%d: %s", aProp->radius,
834  aProp->faccuracy, aProp->niter, aProp->rejsigma,
835  rc, rc != CPL_ERROR_NONE ? cpl_error_get_message() : "");
837 
838  if (!cpl_errorstate_is_equal(prestate)) {
839  cpl_msg_warning(__func__, "The following errors happened while computing "
840  "the astrometric solution (rc = %d/%s):", rc,
841  cpl_error_get_message());
842  cpl_errorstate_dump(prestate, CPL_FALSE, muse_cplerrorstate_dump_some);
843  }
844  return wcs;
845  } /* if MUSE_POSTPROC_ASTROMETRY */
846 
847  rc = muse_wcs_project_tan(pt, aProp->wcs);
848  cpl_msg_debug(__func__, "Astrometry correction returned rc=%d: %s", rc,
849  rc != CPL_ERROR_NONE ? cpl_error_get_message() : "");
850 
851  return pt;
852 } /* muse_postproc_process_exposure() */
853 
854 /*---------------------------------------------------------------------------*/
872 /*---------------------------------------------------------------------------*/
873 cpl_propertylist *
875 {
876  cpl_ensure(aProcessing, CPL_ERROR_NULL_INPUT, NULL);
877 
878  cpl_frameset *fwcs = muse_frameset_find(aProcessing->inframes,
879  MUSE_TAG_OUTPUT_WCS, 0, CPL_FALSE);
880  if (!fwcs || cpl_frameset_get_size(fwcs) <= 0) {
881  /* not an error, quietly return NULL */
882  cpl_frameset_delete(fwcs);
883  return NULL;
884  }
885  cpl_frame *frame = cpl_frameset_get_position(fwcs, 0);
886  const char *fn = cpl_frame_get_filename(frame);
887  cpl_propertylist *header = NULL;
888  int iext, next = cpl_fits_count_extensions(fn);
889  for (iext = 0; iext <= next; iext++) {
890  header = cpl_propertylist_load(fn, iext);
891  cpl_errorstate es = cpl_errorstate_get();
892  cpl_wcs *wcs = cpl_wcs_new_from_propertylist(header);
893  if (!cpl_errorstate_is_equal(es)) {
894  cpl_errorstate_set(es);
895  }
896  if (!wcs) {
897  cpl_propertylist_delete(header);
898  header = NULL;
899  continue; /* try the next extension */
900  }
901  /* need 2D or 3D WCS to have something useful */
902  int naxis = cpl_wcs_get_image_naxis(wcs);
903  cpl_boolean valid = naxis == 2 || naxis == 3;
904  /* need RA---TAN and DEC--TAN for axes 1 and 2 */
905  const cpl_array *ctypes = cpl_wcs_get_ctype(wcs);
906  if (valid && cpl_array_get_string(ctypes, 0)) {
907  valid = valid && !strncmp(cpl_array_get_string(ctypes, 0), "RA---TAN", 9);
908  }
909  if (valid && cpl_array_get_string(ctypes, 1)) {
910  valid = valid && !strncmp(cpl_array_get_string(ctypes, 1), "DEC--TAN", 9);
911  }
912  /* need AWAV or AWAV-LOG for axis 3 */
913  if (valid && cpl_array_get_string(ctypes, 2)) {
914  const char *ctype3 = cpl_array_get_string(ctypes, 2);
915  valid = valid
916  && (!strncmp(ctype3, "AWAV", 5) || !strncmp(ctype3, "AWAV-LOG", 9) ||
917  !strncmp(ctype3, "WAVE", 5) || !strncmp(ctype3, "WAVE-LOG", 9));
918  }
919  if (valid) { /* check for unsupported PCi_j keywords */
920  /* still needs to be done on the input header not on the derives WCS! */
921  cpl_propertylist *pc = cpl_propertylist_new();
922  cpl_propertylist_copy_property_regexp(pc, header, "^PC[0-9]+_[0-9]+", 0);
923  valid = cpl_propertylist_get_size(pc) == 0;
924  cpl_propertylist_delete(pc);
925  }
926  cpl_wcs_delete(wcs);
927  if (valid) {
928  cpl_msg_debug(__func__, "Apparently valid header %s found in extension %d"
929  " of \"%s\"", MUSE_TAG_OUTPUT_WCS, iext, fn);
930  muse_processing_append_used(aProcessing, frame, CPL_FRAME_GROUP_CALIB, 1);
931  break;
932  }
933  cpl_propertylist_delete(header);
934  header = NULL;
935  } /* for iext (all extensions) */
936  if (!header) {
937  cpl_msg_warning(__func__, "No valid headers for %s found in \"%s\"",
938  MUSE_TAG_OUTPUT_WCS, fn);
939  }
940  cpl_frameset_delete(fwcs);
941  return header;
942 } /* muse_postproc_cube_load_output_wcs() */
943 
944 /*---------------------------------------------------------------------------*/
962 /*---------------------------------------------------------------------------*/
963 cpl_error_code
965  muse_pixtable *aPixtable,
966  muse_cube_type aFormat,
967  muse_resampling_params *aParams,
968  const char *aFilter)
969 {
970  cpl_ensure_code(aProcessing && aPixtable && aParams, CPL_ERROR_NULL_INPUT);
971  cpl_ensure_code(aFormat == MUSE_CUBE_TYPE_EURO3D ||
972  aFormat == MUSE_CUBE_TYPE_FITS ||
973  aFormat == MUSE_CUBE_TYPE_EURO3D_X ||
974  aFormat == MUSE_CUBE_TYPE_FITS_X ||
975  aFormat == MUSE_CUBE_TYPE_SDP, CPL_ERROR_ILLEGAL_INPUT);
976 
977  /* convert cube from air to vacuum wavelengths, if needed */
978  if (aParams->tlambda == MUSE_RESAMPLING_DISP_WAVE ||
980  /* convert using standard air (and Ciddor method) */
982  }
983 
984  /* create cube and cast to generic pointer to save code duplication */
985  void *cube = NULL;
986  muse_pixgrid *grid = NULL;
987  if (aFormat == MUSE_CUBE_TYPE_EURO3D || aFormat == MUSE_CUBE_TYPE_EURO3D_X) {
988  cpl_msg_info(__func__, "Resampling to final cube follows, output format "
989  "\"Euro3D\"");
990  cpl_msg_indent_more();
991  muse_euro3dcube *e3d = muse_resampling_euro3d(aPixtable, aParams);
992  cpl_msg_indent_less();
993  cpl_ensure_code(e3d, cpl_error_get_code());
994  cube = e3d;
995  } else if (aFormat == MUSE_CUBE_TYPE_FITS ||
996  aFormat == MUSE_CUBE_TYPE_FITS_X || aFormat == MUSE_CUBE_TYPE_SDP) {
997  cpl_msg_info(__func__, "Resampling to final cube follows, output format "
998  "\"FITS cube\"");
999  cpl_msg_indent_more();
1000  muse_datacube *fits = muse_resampling_cube(aPixtable, aParams, &grid);
1001  cpl_msg_indent_less();
1002  cpl_ensure_code(fits, cpl_error_get_code());
1003  muse_postproc_qc_fwhm(aProcessing, fits);
1004  cube = fits;
1005  }
1006 
1007  cpl_array *filters = muse_cplarray_new_from_delimited_string(aFilter, ",");
1008  int i, ipos = 0, nfilters = cpl_array_get_size(filters);
1009  for (i = 0; i < nfilters; i++) {
1010  /* try to find and load the filter from a table */
1011  muse_table *filter = muse_table_load_filter(aProcessing,
1012  cpl_array_get_string(filters, i));
1013  if (!filter) {
1014  continue;
1015  }
1016 
1017  /* if we got a filter table, do the collapsing */
1018  muse_image *fov = NULL;
1019  if (aFormat == MUSE_CUBE_TYPE_EURO3D || aFormat == MUSE_CUBE_TYPE_EURO3D_X) {
1020  fov = muse_euro3dcube_collapse(cube, filter);
1021  } else if (aFormat == MUSE_CUBE_TYPE_FITS ||
1022  aFormat == MUSE_CUBE_TYPE_FITS_X || aFormat == MUSE_CUBE_TYPE_SDP) {
1023  if (getenv("MUSE_COLLAPSE_PIXTABLE") &&
1024  atoi(getenv("MUSE_COLLAPSE_PIXTABLE")) > 0) {
1025  /* create collapsed image directly from pixel table/grid */
1026  fov = muse_resampling_collapse_pixgrid(aPixtable, grid,
1027  (muse_datacube *)cube,
1028  filter, aParams);
1029  } else {
1030  fov = muse_datacube_collapse(cube, filter);
1031  }
1032  }
1033  muse_table_delete(filter);
1034 
1035  if (aFormat == MUSE_CUBE_TYPE_EURO3D_X || aFormat == MUSE_CUBE_TYPE_FITS_X) {
1036  if (!((muse_datacube *)cube)->recimages) {
1037  /* create empty structures before first entry */
1038  ((muse_datacube *)cube)->recimages = muse_imagelist_new();
1039  ((muse_datacube *)cube)->recnames = cpl_array_new(0, CPL_TYPE_STRING);
1040  }
1041  /* append image to cube/Euro3D structure */
1042  muse_imagelist_set(((muse_datacube *)cube)->recimages, fov,
1043  muse_imagelist_get_size(((muse_datacube *)cube)->recimages));
1044  cpl_array_set_size(((muse_datacube *)cube)->recnames, ipos+1);
1045  cpl_array_set_string(((muse_datacube *)cube)->recnames, ipos,
1046  cpl_array_get_string(filters, i));
1047  } else {
1048  /* we want NANs instead of DQ in the output images */
1049  muse_image_dq_to_nan(fov);
1050  /* duplicate cube header to be able to use its WCS info as the base */
1051  muse_utils_copy_modified_header(fov->header, fov->header, "OBJECT",
1052  cpl_array_get_string(filters, i));
1053  cpl_propertylist_update_string(fov->header, MUSE_HDR_FILTER,
1054  cpl_array_get_string(filters, i));
1055  if (aFormat == MUSE_CUBE_TYPE_SDP) {
1056  if (muse_idp_properties_update_fov(fov) != CPL_ERROR_NONE) {
1057  cpl_msg_warning(__func__, "Writing IDP keywords to field-of-view "
1058  "image failed!");
1059  }
1060  }
1061  muse_processing_save_image(aProcessing, -1, fov, MUSE_TAG_IMAGE_FOV);
1062  muse_image_delete(fov);
1063  }
1064  ipos++;
1065  } /* for i (all filters) */
1066  cpl_array_delete(filters);
1067  muse_pixgrid_delete(grid);
1068 
1069 #ifdef CREATE_IDP
1070  /* Collect IDP properties from the raw data frames and the data cube *
1071  * header */
1072  muse_idp_properties *properties =
1073  muse_idp_properties_collect(aProcessing, ((muse_datacube *)cube)->header,
1074  MUSE_TAG_CUBE_FINAL);
1075 
1076  /* Write formatted IDP properties to the data cube product header so *
1077  * that they are saved to the product file in the next step */
1078  muse_idp_properties_update((muse_datacube *)cube)->header, properties);
1079 #endif
1080 
1081  cpl_error_code rc = CPL_ERROR_NONE;
1082  if (aFormat == MUSE_CUBE_TYPE_EURO3D || aFormat == MUSE_CUBE_TYPE_EURO3D_X) {
1083  /* save the data and clean up */
1084  rc = muse_processing_save_cube(aProcessing, -1, cube, MUSE_TAG_CUBE_FINAL,
1086  muse_euro3dcube_delete(cube);
1087  } else if (aFormat == MUSE_CUBE_TYPE_FITS || aFormat == MUSE_CUBE_TYPE_FITS_X ||
1088  aFormat == MUSE_CUBE_TYPE_SDP) {
1089  /* convert/remove the DQ cube */
1091  /* save the cube and clean up */
1092  if (aFormat == MUSE_CUBE_TYPE_SDP) {
1093  rc = muse_processing_save_cube(aProcessing, -1, cube, MUSE_TAG_CUBE_FINAL,
1095  } else {
1096  rc = muse_processing_save_cube(aProcessing, -1, cube, MUSE_TAG_CUBE_FINAL,
1098  }
1099  muse_datacube_delete(cube);
1100  }
1101 #ifdef CREATE_IDP
1102  muse_idp_properties_delete(properties);
1103 #endif
1104  return rc;
1105 } /* muse_postproc_cube_resample_and_collapse() */
1106 
1107 /*---------------------------------------------------------------------------*/
1130 /*---------------------------------------------------------------------------*/
1131 cpl_error_code
1133 {
1134  cpl_ensure_code(aProcessing && aCube, CPL_ERROR_NULL_INPUT);
1135 
1136  /* determine which FITS keyword prefix to take */
1137  char prefix[KEYWORD_LENGTH] = "";
1138  if (!strncmp(aProcessing->name, "muse_scipost", 13)) {
1139  strncpy(prefix, QC_POST_PREFIX_SCIPOST, KEYWORD_LENGTH);
1140  } else if (!strncmp(aProcessing->name, "muse_exp_combine", 17)) {
1141  strncpy(prefix, QC_POST_PREFIX_EXPCOMBINE, KEYWORD_LENGTH);
1142  } else if (!strncmp(aProcessing->name, "muse_standard", 14)) {
1143  strncpy(prefix, QC_POST_PREFIX_STANDARD, KEYWORD_LENGTH);
1144  } else if (!strncmp(aProcessing->name, "muse_astrometry", 16)) {
1145  strncpy(prefix, QC_POST_PREFIX_ASTROMETRY, KEYWORD_LENGTH);
1146  } else {
1147  cpl_msg_info(__func__, "Recipe \"%s\" found, not generating QC1 keywords",
1148  aProcessing->name);
1149  return CPL_ERROR_NONE;
1150  }
1151 
1152  /* get image from central plane of cube and use it to do object detection *
1153  * and FWHM estimation */
1154  int plane = cpl_imagelist_get_size(aCube->data) / 2;
1155  cpl_image *cim = cpl_imagelist_get(aCube->data, plane);
1156  double dsigmas[] = {/*50., 30., 10., 8., 6.,*/ 5., 4., 3. };
1157  int ndsigmas = sizeof(dsigmas) / sizeof(double);
1158  cpl_vector *vsigmas = cpl_vector_wrap(ndsigmas, dsigmas);
1159  cpl_size isigma = -1;
1160  cpl_errorstate prestate = cpl_errorstate_get();
1161  cpl_apertures *apertures = cpl_apertures_extract(cim, vsigmas, &isigma);
1162  cpl_vector_unwrap(vsigmas);
1163 
1164  /* add the first QC parameters, no matter what the detection did */
1165  char keyword[KEYWORD_LENGTH];
1166  cpl_boolean loglambda = !strncmp(muse_pfits_get_ctype(aCube->header, 3),
1167  "AWAV-LOG", 9);
1168  double crpix3 = muse_pfits_get_crpix(aCube->header, 3),
1169  cd33 = muse_pfits_get_cd(aCube->header, 3, 3),
1170  crval3 = muse_pfits_get_crval(aCube->header, 3),
1171  lambda = loglambda
1172  ? crval3 * exp((plane + 1. - crpix3) * cd33 / crval3)
1173  : (plane + 1. - crpix3) * cd33 + crval3;
1174  snprintf(keyword, KEYWORD_LENGTH, QC_POST_LAMBDA, prefix);
1175  cpl_propertylist_update_float(aCube->header, keyword, lambda);
1176 
1177  /* now see, if we detected anything */
1178  if (!apertures || !cpl_errorstate_is_equal(prestate)) {
1179  /* set up number and a fake entry (source number zero and negative FWHMs) */
1180  snprintf(keyword, KEYWORD_LENGTH, QC_POST_NDET, prefix);
1181  cpl_propertylist_update_int(aCube->header, keyword, 0); /* zero! */
1182 
1183  snprintf(keyword, KEYWORD_LENGTH, QC_POST_POSX, prefix, 0);
1184  cpl_propertylist_update_float(aCube->header, keyword, -1.);
1185  snprintf(keyword, KEYWORD_LENGTH, QC_POST_POSY, prefix, 0);
1186  cpl_propertylist_update_float(aCube->header, keyword, -1.);
1187  snprintf(keyword, KEYWORD_LENGTH, QC_POST_FWHMX, prefix, 0);
1188  cpl_propertylist_update_float(aCube->header, keyword, -1.);
1189  snprintf(keyword, KEYWORD_LENGTH, QC_POST_FWHMY, prefix, 0);
1190  cpl_propertylist_update_float(aCube->header, keyword, -1.);
1191  /* swallow the cpl_aperture error, which is too cryptic for users */
1192  cpl_errorstate_set(prestate);
1193  cpl_msg_warning(__func__, "No sources found for FWHM measurement down to "
1194  "%.1f sigma limit on plane %d, QC parameters will not "
1195  "contain useful information", dsigmas[ndsigmas-1], plane+1);
1196  return CPL_ERROR_DATA_NOT_FOUND;
1197  } /* if no detections */
1198  int ndet = cpl_apertures_get_size(apertures);
1199  snprintf(keyword, KEYWORD_LENGTH, QC_POST_NDET, prefix);
1200  cpl_propertylist_update_int(aCube->header, keyword, ndet); /* zero! */
1201 
1202  /* get some kind of WCS for conversion of FWHM from pixels to arcsec, *
1203  * by default assume WFM and x=RA and y=DEC */
1204  double cd11 = kMuseSpaxelSizeX_WFM, cd12 = 0., cd21 = 0.,
1205  cd22 = kMuseSpaxelSizeY_WFM;
1206  cpl_errorstate es = cpl_errorstate_get();
1207  cpl_wcs *wcs = cpl_wcs_new_from_propertylist(aCube->header);
1208  if (!cpl_errorstate_is_equal(es)) {
1209  cpl_errorstate_set(es); /* swallow possible errors from creating the WCS */
1210  }
1211  if (!wcs || !strncasecmp(muse_pfits_get_ctype(aCube->header, 1), "PIXEL", 5)) {
1212  if (muse_pfits_get_mode(aCube->header) == MUSE_MODE_NFM_AO_N) { /* NFM scaling */
1213  cd11 = kMuseSpaxelSizeX_NFM;
1214  cd22 = kMuseSpaxelSizeY_NFM;
1215  }
1216  } else {
1217  const cpl_matrix *cd = cpl_wcs_get_cd(wcs);
1218  /* take the absolute and scale by 3600 to get positive arcseconds */
1219  cd11 = fabs(cpl_matrix_get(cd, 0, 0)) * 3600.,
1220  cd12 = fabs(cpl_matrix_get(cd, 0, 1)) * 3600.,
1221  cd21 = fabs(cpl_matrix_get(cd, 1, 0)) * 3600.,
1222  cd22 = fabs(cpl_matrix_get(cd, 1, 1)) * 3600.;
1223  }
1224  cpl_wcs_delete(wcs); /* rely on internal NULL check */
1225 
1226  /* Compute FWHM of all apertures; collect them *
1227  * in an image for statistical computations */
1228  cpl_image *statimage = cpl_image_new(ndet, 2, CPL_TYPE_DOUBLE);
1229 #if 0
1230  printf("index RA_FWHM DEC_FWHM\n");
1231 #endif
1232  int n, nbad = 0;
1233  for (n = 1; n <= ndet; n++) {
1234  cpl_size xcen = lround(cpl_apertures_get_centroid_x(apertures, n)),
1235  ycen = lround(cpl_apertures_get_centroid_y(apertures, n));
1236  double x, y;
1237  cpl_errorstate state = cpl_errorstate_get();
1238  cpl_image_get_fwhm(cim, xcen, ycen, &x, &y);
1239  const char *fwhmcomment = NULL;
1240  if (x < 0 || y < 0 || !cpl_errorstate_is_equal(state)) {
1241  /* something's weird with this source */
1242  x = -1.;
1243  y = -1.;
1244  /* do not count these in the statistics */
1245  cpl_image_reject(statimage, n, 1);
1246  cpl_image_reject(statimage, n, 2);
1247  fwhmcomment = "[arcsec] failure determining FWHM";
1248  nbad++;
1249  /* suppress the error state */
1250  cpl_errorstate_set(state);
1251  } else {
1252  /* linear WCS scaling */
1253  x = cd11 * x + cd12 * y;
1254  y = cd22 * y + cd21 * x;
1255  cpl_image_set(statimage, n, 1, x);
1256  cpl_image_set(statimage, n, 2, y);
1257  }
1258 #if 0
1259  printf("%4d %f %f\n", n, x, y);
1260 #endif
1261 
1262  /* generate and write FITS keywords */
1263  snprintf(keyword, KEYWORD_LENGTH, QC_POST_POSX, prefix, n);
1264  cpl_propertylist_update_float(aCube->header, keyword, xcen);
1265  snprintf(keyword, KEYWORD_LENGTH, QC_POST_POSY, prefix, n);
1266  cpl_propertylist_update_float(aCube->header, keyword, ycen);
1267  snprintf(keyword, KEYWORD_LENGTH, QC_POST_FWHMX, prefix, n);
1268  cpl_propertylist_update_float(aCube->header, keyword, x);
1269  if (fwhmcomment) {
1270  cpl_propertylist_set_comment(aCube->header, keyword, fwhmcomment);
1271  }
1272  snprintf(keyword, KEYWORD_LENGTH, QC_POST_FWHMY, prefix, n);
1273  cpl_propertylist_update_float(aCube->header, keyword, y);
1274  if (fwhmcomment) {
1275  cpl_propertylist_set_comment(aCube->header, keyword, fwhmcomment);
1276  }
1277  } /* for n (aperture number) */
1278 #if 0
1279  fflush(stdout);
1280 #endif
1281  cpl_apertures_delete(apertures);
1282 
1283  /* add simple FWHM statistics as well */
1284  snprintf(keyword, KEYWORD_LENGTH, QC_POST_FWHM_NVALID, prefix);
1285  cpl_propertylist_update_int(aCube->header, keyword, ndet - nbad);
1286  cpl_errorstate state = cpl_errorstate_get();
1287  cpl_stats *s = cpl_stats_new_from_image(statimage,
1288  CPL_STATS_MEDIAN | CPL_STATS_MAD);
1289  cpl_boolean error = !cpl_errorstate_is_equal(state);
1290  cpl_errorstate_set(state); /* swallow the error again */
1291  snprintf(keyword, KEYWORD_LENGTH, QC_POST_FWHM_MEDIAN, prefix);
1292  if (error || ndet < 3) {
1293  cpl_propertylist_update_float(aCube->header, keyword, 0.);
1294  } else {
1295  cpl_propertylist_update_float(aCube->header, keyword, cpl_stats_get_median(s));
1296  }
1297  snprintf(keyword, KEYWORD_LENGTH, QC_POST_FWHM_MAD, prefix);
1298  if (error || ndet < 3) {
1299  cpl_propertylist_update_float(aCube->header, keyword, 0.);
1300  } else {
1301  cpl_propertylist_update_float(aCube->header, keyword, cpl_stats_get_mad(s));
1302  }
1303  cpl_stats_delete(s);
1304  cpl_image_delete(statimage);
1305  cpl_msg_info(__func__, "Computed FWHM QC parameters for %d of %d source%s "
1306  "found down to the %.1f sigma threshold on plane %d",
1307  ndet - nbad, ndet, ndet == 1 ? "" : "s", dsigmas[isigma],
1308  plane + 1);
1309 
1310  return CPL_ERROR_NONE;
1311 } /* muse_postproc_qc_fwhm() */
1312 
cpl_error_code muse_sky_lines_fit(cpl_table *, cpl_table *, cpl_image *, muse_wcs *)
Fit all entries of the pixel table to the master sky.
cpl_table * muse_sky_continuum_create(cpl_table *, cpl_table *, cpl_image *, muse_wcs *, double)
Create a continuum spectrum.
muse_image * muse_resampling_collapse_pixgrid(muse_pixtable *aPixtable, muse_pixgrid *aGrid, muse_datacube *aCube, const muse_table *aFilter, muse_resampling_params *aParams)
Integrate a pixel table / pixel grid along the wavelength direction.
#define MUSE_PIXTABLE_FFDATA
Definition: muse_pixtable.h:61
Structure definition of a MUSE datacube.
Definition: muse_datacube.h:48
muse_rvcorrect_type rvtype
Definition: muse_postproc.h:98
cpl_propertylist * muse_postproc_cube_load_output_wcs(muse_processing *aProcessing)
Find a file with a usable output WCS in the input frameset.
muse_postproc_properties * muse_postproc_properties_new(muse_postproc_type aType)
Create a post-processing properties object.
Definition: muse_postproc.c:77
muse_flux_selection_type select
cpl_error_code muse_flux_compute_qc(muse_flux_object *aFluxObj)
Compute QC parameters, related to on-sky throughput.
Definition: muse_flux.c:2858
void muse_image_delete(muse_image *aImage)
Deallocate memory associated to a muse_image object.
Definition: muse_image.c:85
double muse_pfits_get_cd(const cpl_propertylist *aHeaders, unsigned int aAxisI, unsigned int aAxisJ)
find out the WCS coordinate at the reference point
Definition: muse_pfits.c:446
muse_image * intimage
Definition: muse_flux.h:83
muse_pixtable * muse_pixtable_duplicate(muse_pixtable *aPixtable)
Make a copy of the pixtanle.
double muse_pfits_get_crval(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS coordinate at the reference point
Definition: muse_pfits.c:423
#define MUSE_PIXTABLE_FFLAMBDA
Definition: muse_pixtable.h:60
A structure containing a spatial two-axis WCS.
Definition: muse_wcs.h:105
muse_wcs_centroid_type centroid
const char * name
cpl_table * muse_resampling_spectrum_iterate(muse_pixtable *aPixtable, double aBinwidth, float aLo, float aHi, unsigned char aIter)
Iteratively resample selected pixels of a pixel table into spectrum.
double muse_pfits_get_ra(const cpl_propertylist *aHeaders)
find out the right ascension
Definition: muse_pfits.c:267
cpl_image * data
the data extension
Definition: muse_image.h:46
cpl_error_code muse_rvcorrect(muse_pixtable *aPixtable, muse_rvcorrect_type aType)
Correct the wavelengths of all pixels of a given pixel table for radial velocity shift.
void muse_mask_delete(muse_mask *aMask)
Deallocate memory associated to a muse_mask object.
Definition: muse_mask.c:78
muse_xcombine_types muse_postproc_get_weight_type(const char *aWeightString)
Select correct weighting type for weight string.
cpl_error_code muse_dar_correct(muse_pixtable *aPixtable, double aLambdaRef)
Correct the pixel coordinates of all pixels of a given pixel table for differential atmospheric refra...
Definition: muse_dar.c:139
cpl_error_code muse_postproc_revert_ffspec_maybe(muse_pixtable *aPt, const muse_table *aResponse)
Revert correction of on-sky data with the flat-field spectrum.
cpl_image * muse_lsf_average_cube_all(muse_lsf_cube **aLsfCube, muse_pixtable *aPixtable)
Create an average image from all LSF cubes.
cpl_error_code muse_pixtable_and_selected_mask(muse_pixtable *aPixtable, muse_mask *aMask)
Select all pixels where the (x,y) positions are enabled in the given mask.
muse_resampling_crstats_type muse_postproc_get_cr_type(const char *aCRTypeString)
Select correct cosmic ray rejection type for crtype string.
muse_datacube * muse_resampling_cube(muse_pixtable *aPixtable, muse_resampling_params *aParams, muse_pixgrid **aGrid)
Resample a pixel table onto a regular grid structure representing a FITS NAXIS=3 datacube.
void muse_datacube_delete(muse_datacube *aCube)
Deallocate memory associated to a muse_datacube object.
cpl_error_code muse_dar_check(muse_pixtable *aPixtable, double *aMaxShift, cpl_boolean aCorrect, muse_datacube **aCube)
check that the continuum of objects in the frame is well-aligned perpendicular to the spatial axis...
Definition: muse_dar.c:472
#define MUSE_HDR_PT_FFCORR
cpl_error_code muse_flux_reference_table_check(cpl_table *aTable)
Check and/or adapt the standard flux reference table format.
Definition: muse_flux.c:170
cpl_error_code muse_postproc_qc_fwhm(muse_processing *aProcessing, muse_datacube *aCube)
Compute QC1 parameters for datacubes and save them in the FITS header.
void * muse_postproc_process_exposure(muse_postproc_properties *aProp, unsigned int aIndex, muse_postproc_sky_outputs *aSkyOut)
Merge and process pixel tables from one exposure.
WCS object to store data needed while computing the astrometric solution.
Definition: muse_wcs.h:70
The pixel grid.
Definition: muse_pixgrid.h:65
Structure definition of MUSE three extension FITS file.
Definition: muse_image.h:40
muse_flux_object * muse_flux_object_new(void)
Allocate memory for a new muse_flux_object object.
Definition: muse_flux.c:68
cpl_error_code muse_flux_calibrate(muse_pixtable *aPixtable, const muse_table *aResponse, const cpl_table *aExtinction, const muse_table *aTelluric)
Convert the input pixel table from counts to fluxes.
Definition: muse_flux.c:3083
cpl_error_code muse_postproc_cube_resample_and_collapse(muse_processing *aProcessing, muse_pixtable *aPixtable, muse_cube_type aFormat, muse_resampling_params *aParams, const char *aFilter)
High level function to resample to a datacube and collapse that to an image of the field of view and ...
cpl_table * table
The pixel table.
muse_flux_profile_type profile
void muse_lsf_params_delete_all(muse_lsf_params **aParams)
Delete an allocated array of muse_lsf_params structure.
muse_postproc_skymethod skymethod
cpl_propertylist * header
the FITS header
Definition: muse_image.h:72
muse_lsf_params ** lsf_params
unsigned int muse_imagelist_get_size(muse_imagelist *aList)
Return the number of stored images.
cpl_error_code muse_wcs_project_tan(muse_pixtable *aPixtable, const cpl_propertylist *aWCS)
Carry out a gnomonic projection of a pixel table into native spherical coordinates.
Definition: muse_wcs.c:1326
muse_resampling_crstats_type crtype
double muse_pfits_get_crpix(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS reference point
Definition: muse_pfits.c:401
cpl_error_code muse_flux_get_telluric_table(muse_flux_object *aFluxObj)
Get the table of the telluric correction.
Definition: muse_flux.c:2752
void muse_pixgrid_delete(muse_pixgrid *aGrid)
Delete a pixel grid and remove its memory.
Definition: muse_pixgrid.c:575
cpl_array * muse_cplarray_new_from_delimited_string(const char *aString, const char *aDelim)
Convert a delimited string into an array of strings.
Structure definition of MUSE pixel table.
Flux object to store data needed while computing the flux calibration.
Definition: muse_flux.h:78
cpl_error_code muse_datacube_convert_dq(muse_datacube *aCube)
Convert the DQ extension of a datacube to NANs in DATA and STAT.
muse_cube_type muse_postproc_get_cube_format(const char *aFormatString)
Select correct cube format for format string.
cpl_error_code muse_sky_subtract_lines_old(muse_pixtable *aPixtable, cpl_table *aLines, muse_lsf_params **aLsfParams)
Subtract sky lines from a pixtable.
Definition: muse_sky_old.c:532
cpl_array * muse_cpltable_extract_column(cpl_table *aTable, const char *aColumn)
Create an array from a section of a column.
cpl_boolean muse_postproc_check_save_param(const char *aSave, const char *aValid)
Check the –save parameter contents against allowed options.
cpl_error_code muse_flux_response_compute(muse_flux_object *aFluxObj, muse_flux_selection_type aSelect, double aAirmass, const cpl_table *aReference, const cpl_table *aTellBands, const cpl_table *aExtinct)
Compare measured flux distribution over wavelength with calibrated stellar fluxes and derive instrume...
Definition: muse_flux.c:2070
muse_resampling_crstats_type
Cosmic ray rejection statistics type.
muse_resampling_params * muse_resampling_params_new(muse_resampling_type aMethod)
Create the resampling parameters structure.
cpl_error_code muse_processing_save_cube(muse_processing *aProcessing, int aIFU, void *aCube, const char *aTag, muse_cube_type aType)
Save a MUSE datacube to disk.
cpl_error_code muse_wcs_solve(muse_wcs_object *aWCSObj, cpl_table *aReference, float aRadius, float aFAccuracy, int aIter, float aThresh)
Find the world coordinate solution for a given field using coordinates of detected sources and coordi...
Definition: muse_wcs.c:671
muse_wcs_object * muse_wcs_object_new(void)
Allocate memory for a new muse_wcs_object object.
Definition: muse_wcs.c:71
cpl_table * ffspec
A flat-field spectrum.
muse_flux_smooth_type smooth
muse_cube_type
cpl_propertylist * wcs
Structure definition of the post-processing properties.
Definition: muse_postproc.h:91
void muse_cplerrorstate_dump_some(unsigned aCurrent, unsigned aFirst, unsigned aLast)
Dump some CPL errors.
void muse_euro3dcube_delete(muse_euro3dcube *aEuro3D)
Deallocate memory associated to a muse_euro3dcube object.
muse_table * muse_table_load_filter(muse_processing *aProcessing, const char *aFilterName)
Load a table for a given filter name.
Definition: muse_utils.c:799
void muse_postproc_properties_delete(muse_postproc_properties *aProp)
Free memory taken by a post-processing properties object and all its components.
cpl_table * muse_postproc_load_nearest(const cpl_propertylist *aHeader, const cpl_frame *aFrame, float aWarnLimit, float aErrLimit, double *aRA, double *aDEC)
Load the calibration from a multi-table FITS file that is nearest on the sky.
cpl_imagelist * data
the cube containing the actual data values
Definition: muse_datacube.h:76
cpl_error_code muse_sky_subtract_continuum(muse_pixtable *, cpl_table *)
Subtract sky continuum from pixel table.
cpl_error_code muse_sky_lines_set_range(cpl_table *, double, double)
Limit the lines in the table to a wavelength range.
Structure definition of a Euro3D datacube.
Definition: muse_datacube.h:97
double muse_pfits_get_dec(const cpl_propertylist *aHeaders)
find out the declination
Definition: muse_pfits.c:285
void muse_processing_append_used(muse_processing *aProcessing, cpl_frame *aFrame, cpl_frame_group aGroup, int aDuplicate)
Add a frame to the set of used frames.
Structure to store a table together with a property list.
Definition: muse_table.h:43
cpl_error_code muse_flux_get_response_table(muse_flux_object *aFluxObj, muse_flux_smooth_type aSmooth)
Get the table of the standard star response function.
Definition: muse_flux.c:2644
void muse_table_delete(muse_table *aTable)
Deallocate memory associated to a muse_table object.
Definition: muse_table.c:80
muse_resampling_dispersion_type tlambda
muse_lsf_cube ** lsf_cube
cpl_propertylist * header
the header
Definition: muse_table.h:56
const char * muse_pfits_get_ctype(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS axis type string
Definition: muse_pfits.c:469
muse_postproc_type
Type of per-exposure processing to run.
Definition: muse_postproc.h:58
cpl_error_code muse_flux_integrate_std(muse_pixtable *aPixtable, muse_flux_profile_type aProfile, muse_flux_object *aFluxObj)
Reconstruct a cube, detect the standard star, and integrate its flux.
Definition: muse_flux.c:1362
cpl_propertylist * header
the FITS header
Definition: muse_datacube.h:57
int muse_processing_save_image(muse_processing *aProcessing, int aIFU, muse_image *aImage, const char *aTag)
Save a computed MUSE image to disk.
muse_pixtable * muse_pixtable_load_merge_channels(cpl_table *aExposureList, double aLambdaMin, double aLambdaMax)
Load and merge the pixel tables of the 24 MUSE sub-fields.
cpl_error_code muse_wcs_optimize_solution(muse_wcs_object *aWCSObj, float aDetSigma, muse_wcs_centroid_type aCentroid, cpl_table *aReference, float aRadius, float aFAccuracy, int aIter, float aRejSigma)
Find an optimal astrometry solution given a detection image.
Definition: muse_wcs.c:988
muse_resampling_type muse_postproc_get_resampling_type(const char *aResampleString)
Select correct resampling type for resample string.
cpl_error_code muse_phys_air_to_vacuum(muse_pixtable *aPixtable, unsigned int aFlags)
Convert a pixel table from air to vacuum wavelengths.
Definition: muse_phys.c:261
Handling of "mask" files.
Definition: muse_mask.h:43
void muse_lsf_cube_delete_all(muse_lsf_cube **aLsfCube)
Delete all LSF cubes.
muse_imagelist * muse_imagelist_new(void)
Create a new (empty) MUSE image list.
cpl_frameset * inframes
muse_mask * muse_sky_create_skymask(muse_image *, double, double, const char *)
Select spaxels to be considered as sky.
cpl_error_code muse_sky_subtract_lines(muse_pixtable *, cpl_table *, muse_lsf_cube **)
Subtract sky lines from a pixtable.
muse_postproc_type type
Definition: muse_postproc.h:92
muse_resampling_type
Resampling types.
cpl_error_code muse_pixtable_spectrum_apply(muse_pixtable *aPixtable, const cpl_array *aLambdas, const cpl_array *aFlux, muse_pixtable_operation aOperation)
Apply a spectrum given by two arrays with an operation to a pixel table.
cpl_error_code muse_lsf_fold_rectangle(cpl_image *aLsfImage, const muse_wcs *aWCS, double aBinWidth)
Filter an LSF image with a rectangle to model spectrum binning.
muse_euro3dcube * muse_resampling_euro3d(muse_pixtable *aPixtable, muse_resampling_params *aParams)
Resample a pixel table onto a regular grid structure representing a Euro3D format file...
Resampling parameters.
Structure definition of the post-processing output sky data.
void muse_resampling_params_delete(muse_resampling_params *aParams)
Delete a resampling parameters structure.
muse_image * muse_datacube_collapse(muse_datacube *aCube, const muse_table *aFilter)
Integrate a FITS NAXIS=3 datacube along the wavelength direction.
cpl_table * muse_resampling_spectrum(muse_pixtable *aPixtable, double aBinwidth)
Resample the selected pixels of a pixel table into a spectrum.
cpl_frameset * muse_frameset_find(const cpl_frameset *aFrames, const char *aTag, unsigned char aIFU, cpl_boolean aInvert)
return frameset containing data from an IFU/channel with a certain tag
Definition: muse_utils.c:160
cpl_error_code muse_sky_lines_fit_old(cpl_table *aSpectrum, cpl_table *aLines)
Fit all entries of the pixel table to the master sky.
Definition: muse_sky_old.c:436
muse_sky_params skymodel_params
void muse_pixtable_delete(muse_pixtable *aPixtable)
Deallocate memory associated to a pixel table object.
muse_postproc_darcheck darcheck
Definition: muse_postproc.h:97
muse_xcombine_types
Xposure combination types.
Definition: muse_xcombine.h:44
cpl_error_code muse_imagelist_set(muse_imagelist *aList, muse_image *aImage, unsigned int aIdx)
Set the muse_image of given list index.
double muse_astro_angular_distance(double aRA1, double aDEC1, double aRA2, double aDEC2)
Compute angular distance in the sky between two positions.
Definition: muse_astro.c:496
muse_ins_mode muse_pfits_get_mode(const cpl_propertylist *aHeaders)
find out the observation mode
Definition: muse_pfits.c:1352
cpl_error_code muse_utils_copy_modified_header(cpl_propertylist *aH1, cpl_propertylist *aH2, const char *aKeyword, const char *aString)
Copy a modified header keyword from one header to another.
Definition: muse_utils.c:2570
cpl_error_code muse_wcs_locate_sources(muse_pixtable *aPixtable, float aSigma, muse_wcs_centroid_type aCentroid, muse_wcs_object *aWCSObj)
Determine the centers of stars on an astrometric exposure.
Definition: muse_wcs.c:496
muse_image * muse_euro3dcube_collapse(muse_euro3dcube *aEuro3D, const muse_table *aFilter)
Integrate a Euro3D datacube along the wavelength direction.
cpl_error_code muse_image_dq_to_nan(muse_image *aImage)
Convert pixels flagged in the DQ extension to NANs in DATA (and STAT, if present).
Definition: muse_image.c:904
double muse_astro_airmass(cpl_propertylist *aHeader)
Derive the effective airmass of an observation from information in a FITS header. ...
Definition: muse_astro.c:331
cpl_propertylist * header
The FITS header.