MUSE Pipeline Reference Manual  2.1.1
muse_pixtable.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-2016 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 <math.h>
31 #include <string.h>
32 
33 #include "muse_pixtable.h"
34 #include "muse_instrument.h"
35 
36 #include "muse_cplwrappers.h"
37 #include "muse_geo.h"
38 #include "muse_mask.h"
39 #include "muse_pfits.h"
40 #include "muse_quadrants.h"
41 #include "muse_quality.h"
42 #include "muse_tracing.h"
43 #include "muse_wavecalib.h"
44 #include "muse_wcs.h"
45 #include "muse_utils.h"
46 
47 /*----------------------------------------------------------------------------*
48  * Debugging Macros *
49  * Set these to 1 or higher for (lots of) debugging output *
50  *----------------------------------------------------------------------------*/
51 #define CREATE_MINIMAL_PIXTABLE 0 /* create a very small pixel table */
52 #define PIXTABLE_CREATE_CCDSIZED 1 /* start with pixel table sized like input image */
53 #define DEBUG_PIXTABLE_CREATION 0 /* some table dumps in muse_pixtable_create() */
54 #define DEBUG_PIXTABLE_FEW_SLICES 0 /* when doing quick pixel table creation for *
55  * debugging, set this to a value below 48 */
56 
57 /*----------------------------------------------------------------------------*/
95 /*----------------------------------------------------------------------------*/
96 
99 /* Let's use macros for the shift lengths we need */
100 #define MUSE_ORIGIN_SHIFT_XSLICE 24
101 #define MUSE_ORIGIN_SHIFT_YPIX 11
102 #define MUSE_ORIGIN_SHIFT_IFU 6
103 /* spacing between the left slice edge at the bottom of the slice *
104  * to the origin of the slice from where we could x positions within *
105  * the slice; this is to make sure that no pixel is left of the offset */
106 #define MUSE_ORIGIN_SLICE_SAFETY_OFFSET -20
107 
108 /* private functions without error checking */
110 static inline uint32_t
111 muse_pixtable_origin_encode_fast(unsigned int aX, unsigned int aY,
112  unsigned short aIFU,
113  unsigned short aSlice, unsigned int aOffset)
114 {
115  return ((aX - aOffset) << MUSE_ORIGIN_SHIFT_XSLICE)
116  | (aY << MUSE_ORIGIN_SHIFT_YPIX)
117  | (aIFU << MUSE_ORIGIN_SHIFT_IFU)
118  | aSlice;
119 }
120 
122 static inline unsigned int
123 muse_pixtable_origin_get_x_fast(uint32_t aOrigin, uint32_t aOffset)
124 {
125  return ((aOrigin >> MUSE_ORIGIN_SHIFT_XSLICE) & 0x7f) + aOffset;
126 }
127 
129 static inline unsigned int
130 muse_pixtable_origin_get_y_fast(uint32_t aOrigin)
131 {
132  return (aOrigin >> MUSE_ORIGIN_SHIFT_YPIX) & 0x1fff;
133 }
134 
136 static inline unsigned short
137 muse_pixtable_origin_get_ifu_fast(uint32_t aOrigin)
138 {
139  return (aOrigin >> MUSE_ORIGIN_SHIFT_IFU) & 0x1f;
140 }
141 
143 static inline unsigned short
144 muse_pixtable_origin_get_slice_fast(uint32_t aOrigin)
145 {
146  return aOrigin & 0x3f;
147 }
148 
149 /*---------------------------------------------------------------------------*/
165 /*---------------------------------------------------------------------------*/
166 uint32_t
167 muse_pixtable_origin_encode(unsigned int aX, unsigned int aY,
168  unsigned short aIFU,
169  unsigned short aSlice, unsigned int aOffset)
170  /* use the explicit 32bit return type so that *
171  * this also works on systems with other setups */
172 {
173  /* check for allowed values to fit into 32bit */
174  cpl_ensure(aX < 8192 && aX > 0 && aY < 8192 && aY > 0 &&
175  aIFU <= kMuseNumIFUs && aIFU >= 1 &&
176  aSlice <= kMuseSlicesPerCCD && aSlice >= 1 && aOffset < 8192,
177  CPL_ERROR_ILLEGAL_INPUT, 0);
178 
179 #if 0
180  cpl_msg_debug(__func__, "origin (%d, %d, %d, %d, %d) = 0x%x",
181  aX, aY, aIFU, aSlice, aOffset,
182  ((aX - aOffset) << MUSE_ORIGIN_SHIFT_XSLICE)
183  | (aY << MUSE_ORIGIN_SHIFT_YPIX)
184  | (aIFU << MUSE_ORIGIN_SHIFT_IFU)
185  | aSlice);
186 #endif
187 
188  /* do the bit-shifting "magic" */
189  return muse_pixtable_origin_encode_fast(aX, aY, aIFU, aSlice, aOffset);
190 } /* muse_pixtable_origin_encode() */
191 
192 /*---------------------------------------------------------------------------*/
210 /*---------------------------------------------------------------------------*/
211 unsigned int
212 muse_pixtable_origin_get_x(uint32_t aOrigin, muse_pixtable *aPixtable,
213  cpl_size aRow)
214 {
215  unsigned short slice = muse_pixtable_origin_get_slice_fast(aOrigin),
216  ifu = muse_pixtable_origin_get_ifu_fast(aOrigin);
217  cpl_errorstate prestate = cpl_errorstate_get();
218  unsigned int expnum = muse_pixtable_get_expnum(aPixtable, aRow);
219  if (!cpl_errorstate_is_equal(prestate)) {
220  cpl_errorstate_set(prestate);
221  }
222  unsigned int offset = muse_pixtable_origin_get_offset(aPixtable, expnum, ifu,
223  slice),
224  x = muse_pixtable_origin_get_x_fast(aOrigin, offset);
225 #if 0
226  if (x > 8191 || x < 1 || !cpl_errorstate_is_equal(prestate)) {
227  cpl_msg_error(__func__, "aOrigin=%#x x=%d (%d %d %d), %s",
228  aOrigin, x, slice, ifu, offset, cpl_error_get_message());
229  }
230 #endif
231  cpl_ensure(x <= 8191 && x >= 1 && cpl_errorstate_is_equal(prestate),
232  CPL_ERROR_ILLEGAL_OUTPUT, 0);
233  return x;
234 } /* muse_pixtable_origin_get_x() */
235 
236 /*---------------------------------------------------------------------------*/
244 /*---------------------------------------------------------------------------*/
245 unsigned int
246 muse_pixtable_origin_get_y(uint32_t aOrigin)
247 {
248  unsigned int y = muse_pixtable_origin_get_y_fast(aOrigin);
249 #if 0
250  if (y > 8191 || y < 1) {
251  cpl_msg_error(__func__, "aOrigin=%#x y=%d", aOrigin, y);
252  }
253 #endif
254  cpl_ensure(y <= 8191 && y >= 1, CPL_ERROR_ILLEGAL_OUTPUT, 0);
255  return y;
256 } /* muse_pixtable_origin_get_y() */
257 
258 /*---------------------------------------------------------------------------*/
266 /*---------------------------------------------------------------------------*/
267 unsigned short
269 {
270  unsigned short ifu = muse_pixtable_origin_get_ifu_fast(aOrigin);
271 #if 0
272  if (ifu > kMuseNumIFUs || ifu < 1) {
273  cpl_msg_error(__func__, "aOrigin=%#x ifu=%d", aOrigin, ifu);
274  }
275 #endif
276  cpl_ensure(ifu <= kMuseNumIFUs && ifu >= 1, CPL_ERROR_ILLEGAL_OUTPUT, 0);
277  return ifu;
278 } /* muse_pixtable_origin_get_ifu() */
279 
280 /*---------------------------------------------------------------------------*/
288 /*---------------------------------------------------------------------------*/
289 unsigned short
291 {
292  unsigned short slice = muse_pixtable_origin_get_slice_fast(aOrigin);
293 #if 0
294  if (slice > kMuseSlicesPerCCD || slice < 1) {
295  cpl_msg_error(__func__, "aOrigin=%#x slice=%d", aOrigin, slice);
296  }
297 #endif
298  cpl_ensure(slice <= kMuseSlicesPerCCD && slice >= 1,
299  CPL_ERROR_ILLEGAL_OUTPUT, 0);
300  return slice;
301 } /* muse_pixtable_origin_get_slice() */
302 
303 /*---------------------------------------------------------------------------*/
318 /*---------------------------------------------------------------------------*/
319 unsigned int
321  cpl_polynomial *aLTrace,
322  unsigned short aIFU, unsigned short aSlice)
323 {
324  cpl_ensure(aPixtable && aPixtable->header, CPL_ERROR_NULL_INPUT, 0);
325  cpl_errorstate prestate = cpl_errorstate_get();
326  unsigned int offset = floor(cpl_polynomial_eval_1d(aLTrace, 1, NULL))
327  + MUSE_ORIGIN_SLICE_SAFETY_OFFSET;
328  cpl_ensure(cpl_errorstate_is_equal(prestate), cpl_error_get_code(), 0);
329  /* use exposure number 0 here, this is set only later, *
330  * when merging pixel tables of different exposures */
331  char *keyword = cpl_sprintf(MUSE_HDR_PT_IFU_SLICE_OFFSET, 0, aIFU, aSlice);
332  cpl_propertylist_update_int(aPixtable->header, keyword, offset);
333  cpl_propertylist_set_comment(aPixtable->header, keyword,
334  MUSE_HDR_PT_IFU_SLICE_OFFSET_COMMENT);
335  cpl_free(keyword);
336  return offset;
337 } /* muse_pixtable_origin_set_offset() */
338 
339 /*---------------------------------------------------------------------------*/
353 /*---------------------------------------------------------------------------*/
354 unsigned int
355 muse_pixtable_origin_get_offset(muse_pixtable *aPixtable, unsigned int aExpNum,
356  unsigned short aIFU, unsigned short aSlice)
357 {
358  cpl_ensure(aPixtable && aPixtable->header, CPL_ERROR_NULL_INPUT, 0);
359  char *keyword = cpl_sprintf(MUSE_HDR_PT_IFU_SLICE_OFFSET, aExpNum, aIFU,
360  aSlice);
361  cpl_errorstate prestate = cpl_errorstate_get();
362  unsigned int offset = cpl_propertylist_get_int(aPixtable->header, keyword);
363  cpl_free(keyword);
364 #if 0
365  if (offset > 8191 || offset < 1) {
366  cpl_msg_error(__func__, "aIFU=%d aSlice=%d offset=%d",
367  aIFU, aSlice, offset);
368  }
369 #endif
370  cpl_ensure(offset <= 8191 && offset >= 1 && cpl_errorstate_is_equal(prestate),
371  CPL_ERROR_ILLEGAL_OUTPUT, 0);
372  return offset;
373 } /* muse_pixtable_origin_get_offset() */
374 
375 /*---------------------------------------------------------------------------*/
395 /*---------------------------------------------------------------------------*/
396 cpl_error_code
398  unsigned int aNum)
399 {
400  cpl_ensure_code(aOut && aOut->header, CPL_ERROR_NULL_INPUT);
401  cpl_propertylist *dest = aOut->header,
402  *from = aOut->header;
403  if (aFrom && aFrom->header) {
404  from = aFrom->header;
405  }
406  char keyword[KEYWORD_LENGTH];
407  unsigned short nifu, nslice;
408  for (nifu = 1; nifu <= kMuseNumIFUs; nifu++) {
409  for (nslice = 1; nslice <= kMuseSlicesPerCCD; nslice++) {
410  /* keyword to copy */
411  snprintf(keyword, KEYWORD_LENGTH, MUSE_HDR_PT_IFU_SLICE_OFFSET,
412  0, nifu, nslice);
413  cpl_errorstate prestate = cpl_errorstate_get();
414  unsigned int offset = cpl_propertylist_get_int(from, keyword);
415  if (!cpl_errorstate_is_equal(prestate)) {
416  /* this one was not set, skip to the next one */
417  cpl_errorstate_set(prestate);
418  continue;
419  }
420  /* erase old keyword */
421  cpl_propertylist_erase(from, keyword);
422  /* add new keyword */
423  snprintf(keyword, KEYWORD_LENGTH, MUSE_HDR_PT_IFU_SLICE_OFFSET,
424  aNum, nifu, nslice);
425  cpl_propertylist_update_int(dest, keyword, offset);
426  cpl_propertylist_set_comment(dest, keyword,
427  MUSE_HDR_PT_IFU_SLICE_OFFSET_COMMENT);
428  } /* for nslice */
429  } /* for nifu */
430 
431  return CPL_ERROR_NONE;
432 } /* muse_pixtable_origin_copy_offsets() */
433 
434 /*---------------------------------------------------------------------------*/
449 /*---------------------------------------------------------------------------*/
450 unsigned int
451 muse_pixtable_get_expnum(muse_pixtable *aPixtable, cpl_size aRow)
452 {
453  cpl_ensure(aPixtable && aPixtable->header, CPL_ERROR_NULL_INPUT, 0);
454  cpl_ensure(aRow >= 0 && muse_pixtable_get_nrow(aPixtable) > aRow,
455  CPL_ERROR_ILLEGAL_INPUT, 0);
456 
457  char keyword[KEYWORD_LENGTH];
458  unsigned int exposure = 0;
459  cpl_size lo = 0, hi = 0;
460  do {
461  cpl_errorstate prestate = cpl_errorstate_get();
462  snprintf(keyword, KEYWORD_LENGTH, MUSE_HDR_PT_EXP_FST, ++exposure);
463  lo = cpl_propertylist_get_long_long(aPixtable->header, keyword);
464  snprintf(keyword, KEYWORD_LENGTH, MUSE_HDR_PT_EXP_LST, exposure);
465  hi = cpl_propertylist_get_long_long(aPixtable->header, keyword);
466  if (!cpl_errorstate_is_equal(prestate)) {
467  if (exposure == 1) {
468  /* no exposure headers in the table at all, so it's either the *
469  * only one or the table was sorted, destroying exposure info */
470  lo = hi = aRow;
471  exposure = 0; /* signify the problem to the caller */
472  }
473  cpl_errorstate_set(prestate);
474  break;
475  }
476  } while (hi < aRow);
477  cpl_ensure(lo <= aRow && hi >= aRow, CPL_ERROR_ILLEGAL_OUTPUT, 0);
478  return exposure;
479 } /* muse_pixtable_get_expnum */
480 
481 /*----------------------------------------------------------------------------*/
501 /*----------------------------------------------------------------------------*/
503  { MUSE_PIXTABLE_XPOS, CPL_TYPE_FLOAT, "pix", "%7.2f",
504  "relative x-pixel position in the output datacube", CPL_TRUE},
505  { MUSE_PIXTABLE_YPOS, CPL_TYPE_FLOAT, "pix", "%7.2f",
506  "relative y-pixel position in the output datacube", CPL_TRUE},
507  { MUSE_PIXTABLE_LAMBDA, CPL_TYPE_FLOAT, "Angstrom", "%8.2f",
508  "wavelength of this pixel", CPL_TRUE},
509  { MUSE_PIXTABLE_DATA, CPL_TYPE_FLOAT, "count", "%e",
510  "data value in this pixel", CPL_TRUE},
511  { MUSE_PIXTABLE_DQ, CPL_TYPE_INT, NULL, "%#x",
512  "Euro3D bad pixel status of this pixel", CPL_TRUE},
513  { MUSE_PIXTABLE_STAT, CPL_TYPE_FLOAT, "count**2", "%e",
514  "data variance of this pixel", CPL_TRUE},
515  { MUSE_PIXTABLE_ORIGIN, CPL_TYPE_INT, NULL, "0x%08x",
516  "encoded value of IFU and slice number, as well as x and y "
517  "position in the raw (trimmed) data", CPL_TRUE},
518  { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
519 };
520 
521 /*---------------------------------------------------------------------------*/
567 /*---------------------------------------------------------------------------*/
569 muse_pixtable_create(muse_image *aImage, cpl_table *aTrace, cpl_table *aWave,
570  cpl_table *aGeoTable)
571 {
572  cpl_ensure(aImage && aWave && aTrace, CPL_ERROR_NULL_INPUT, NULL);
573  /* get and check IFU number */
574  unsigned char ifu = muse_utils_get_ifu(aImage->header);
575  cpl_ensure(ifu >= 1 && ifu <= kMuseNumIFUs, CPL_ERROR_DATA_NOT_FOUND, NULL);
576 
577  /* check that the input data is in count units */
578  if (!cpl_propertylist_has(aImage->header, "BUNIT")) {
579  char *msg = cpl_sprintf("Input data of IFU %hhu does not contain a data "
580  "unit (\"BUNIT\"), cannot proceed!", ifu);
581  cpl_msg_error(__func__, "%s", msg);
582  cpl_error_set_message(__func__, CPL_ERROR_INCOMPATIBLE_INPUT, "%s", msg);
583  cpl_free(msg);
584  return NULL;;
585  }
586  const char *unit = muse_pfits_get_bunit(aImage->header);
587  if (strncmp(unit, "count", 6)) {
588  char *msg = cpl_sprintf("Input data of IFU %hhu is not in \"count\" units "
589  "but in \"%s\", not suitable for a new pixel table!",
590  ifu, unit);
591  cpl_msg_error(__func__, "%s", msg);
592  cpl_error_set_message(__func__, CPL_ERROR_INCOMPATIBLE_INPUT, "%s", msg);
593  cpl_free(msg);
594  return NULL;
595  }
596 
597  /* create the output pixel table structure, and NULL out the components */
598  muse_pixtable *pt = cpl_calloc(1, sizeof(muse_pixtable));
599  /* compute initial table size */
600  int xsize = cpl_image_get_size_x(aImage->data),
601  ysize = cpl_image_get_size_y(aImage->data),
602 #if PIXTABLE_CREATE_CCDSIZED
603  /* There cannot be more pixels needed in the table than in the *
604  * input image, so with this we can remove the overflow check below. */
605  isize = xsize * ysize;
606 #else
607  /* It's very likely that no more pixels than these are needed. But then *
608  * we need to check for overflows below (and still cut it in the end). */
609  isize = ysize * kMuseSlicesPerCCD * kMuseSliceHiLikelyWidth;
610 #endif
611  pt->table = muse_cpltable_new(muse_pixtable_def, isize);
612 
613  /* copy the header from the input image to the pixel table, *
614  * but remove some keywords that could cause trouble */
615  pt->header = cpl_propertylist_duplicate(aImage->header);
616  cpl_propertylist_erase_regexp(pt->header,
617  "^SIMPLE$|^BITPIX$|^NAXIS|^EXTEND$|^XTENSION$|"
618  "^DATASUM$|^DATAMIN$|^DATAMAX$|^DATAMD5$|"
619  "^PCOUNT$|^GCOUNT$|^HDUVERS$|^BLANK$|"
620  "^BZERO$|^BSCALE$|^BUNIT$|^CHECKSUM$|^INHERIT$|"
621  MUSE_HDR_OVSC_REGEXP"|"MUSE_WCS_KEYS, 0);
622  cpl_propertylist_append_string(pt->header, MUSE_HDR_PT_TYPE,
623  aGeoTable ? MUSE_PIXTABLE_STRING_FULL
624  : MUSE_PIXTABLE_STRING_SIMPLE);
625  cpl_propertylist_set_comment(pt->header, MUSE_HDR_PT_TYPE,
626  aGeoTable ? MUSE_PIXTABLE_COMMENT_FULL
627  : MUSE_PIXTABLE_COMMENT_SIMPLE);
628 
629  /* make all table cells valid by filling them with zeros; this is a workaround *
630  * suggested by Carlo Izzo to be able to set the table cells using array *
631  * access instead of the slower and non-threadsafe cpl_table_set() function */
632  cpl_table_fill_column_window_float(pt->table, MUSE_PIXTABLE_XPOS, 0, isize, 0);
633  cpl_table_fill_column_window_float(pt->table, MUSE_PIXTABLE_YPOS, 0, isize, 0);
634  cpl_table_fill_column_window_float(pt->table, MUSE_PIXTABLE_LAMBDA, 0, isize, 0);
635  cpl_table_fill_column_window_float(pt->table, MUSE_PIXTABLE_DATA, 0, isize, 0);
636  cpl_table_fill_column_window_float(pt->table, MUSE_PIXTABLE_STAT, 0, isize, 0);
637  cpl_table_fill_column_window_int(pt->table, MUSE_PIXTABLE_DQ, 0, isize, 0);
638  cpl_table_fill_column_window_int(pt->table, MUSE_PIXTABLE_ORIGIN, 0, isize, 0);
639 
640  /* get columns as pointers */
641  float *cdata_xpos = cpl_table_get_data_float(pt->table, MUSE_PIXTABLE_XPOS),
642  *cdata_ypos = cpl_table_get_data_float(pt->table, MUSE_PIXTABLE_YPOS),
643  *cdata_lambda = cpl_table_get_data_float(pt->table, MUSE_PIXTABLE_LAMBDA),
644  *cdata_data = cpl_table_get_data_float(pt->table, MUSE_PIXTABLE_DATA),
645  *cdata_stat = cpl_table_get_data_float(pt->table, MUSE_PIXTABLE_STAT);
646  int *cdata_dq = cpl_table_get_data_int(pt->table, MUSE_PIXTABLE_DQ);
647  uint32_t *cdata_origin = (uint32_t *)cpl_table_get_data_int(pt->table,
649 
650  /* get input image buffers as pointers, too */
651  const float *pixdata = cpl_image_get_data_float_const(aImage->data),
652  *pixstat = cpl_image_get_data_float_const(aImage->stat);
653  const int *pixdq = cpl_image_get_data_int_const(aImage->dq);
654 
655  /* get columns of the wavecal table */
656  unsigned short wavexorder, waveyorder;
657  muse_wave_table_get_orders(aWave, &wavexorder, &waveyorder);
658  cpl_msg_info(__func__, "Creating pixel table for IFU %hhu, using order %d for "
659  "trace solution and orders %hu/%hu for wavelength solution", ifu,
660  muse_trace_table_get_order(aTrace), wavexorder, waveyorder);
661 
662  /* get position corresponding to the IFU number */
663  cpl_table *geopos = NULL;
664  double *slice_x = NULL, *slice_y = NULL, *slice_angle = NULL,
665  *slice_width = NULL;
666  if (aGeoTable) {
667  geopos = muse_geo_table_extract_ifu(aGeoTable, ifu);
668  slice_x = cpl_table_get_data_double(geopos, MUSE_GEOTABLE_X);
669  slice_y = cpl_table_get_data_double(geopos, MUSE_GEOTABLE_Y);
670  /* fail on missing x/y as they are critical */
671  if (!geopos || !slice_x || !slice_y) {
672  char *msg = !geopos
673  ? cpl_sprintf("Geometry table is missing data for IFU %hhu!", ifu)
674  : cpl_sprintf("Geometry table is missing column%s%s%s!",
675  !slice_x && !slice_y ? "s" : "",
676  slice_x ? "" : " \""MUSE_GEOTABLE_X"\"",
677  slice_y ? "" : " \""MUSE_GEOTABLE_Y"\"");
678  cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT, "%s", msg);
679  cpl_msg_error(__func__, "%s", msg);
680  cpl_free(msg);
681  cpl_table_delete(geopos);
683  return NULL;
684  }
685  /* only output warnings for other missing columns, they *
686  * default to usable values below, ignore CPL errors */
687  cpl_errorstate prestate = cpl_errorstate_get();
688  slice_angle = cpl_table_get_data_double(geopos, MUSE_GEOTABLE_ANGLE);
689  slice_width = cpl_table_get_data_double(geopos, MUSE_GEOTABLE_WIDTH);
690  if (!cpl_errorstate_is_equal(prestate)) {
691  char *msg = cpl_sprintf("Geometry table is missing column%s%s%s!",
692  !slice_angle && !slice_width ? "s" : "",
693  slice_angle ? "" : " \""MUSE_GEOTABLE_ANGLE"\"",
694  slice_width ? "" : " \""MUSE_GEOTABLE_WIDTH"\"");
695  cpl_error_set_message(__func__, CPL_ERROR_BAD_FILE_FORMAT, "%s", msg);
696  cpl_msg_error(__func__, "%s", msg);
697  cpl_free(msg);
698  if (!getenv("MUSE_EXPERT_USER")) {
699  cpl_table_delete(geopos);
701  return NULL;
702  }
703  } /* if */
704  } /* if aGeoTable */
705 
706  /* this counts the current table row that is being entered */
707  cpl_size itablerow = 0;
708  /* time the operation for this IFU */
709  double cputime = cpl_test_get_cputime(),
710  walltime = cpl_test_get_walltime();
711  /* loop through all slices */
712  int islice;
713 #if DEBUG_PIXTABLE_FEW_SLICES > 0
714 # define SLICE_LIMIT DEBUG_PIXTABLE_FEW_SLICES
715 #else
716 # define SLICE_LIMIT kMuseSlicesPerCCD
717 #endif
718  for (islice = 0; islice < SLICE_LIMIT; islice++) {
719 #if DEBUG_PIXTABLE_CREATION /* check how many pixels extend over the nominal width */
720  cpl_msg_debug(__func__, "Starting to process slice %2d of IFU %2hhu",
721  islice + 1, ifu);
722 #endif
723 
724  /* fill the wavelength calibration polynomial for this slice */
725  cpl_polynomial *pwave = muse_wave_table_get_poly_for_slice(aWave, islice + 1);
726  /* vector for the position within the slice (for evaluation *
727  * of the wavelength solution in two dimensions */
728  cpl_vector *pos = cpl_vector_new(2);
729 
730  /* test evaluate the wavelength polynomial to see if it's valid */
731  cpl_vector_set(pos, 0, kMuseOutputYTop / 2);
732  cpl_vector_set(pos, 1, kMuseOutputYTop / 2);
733  double ltest = cpl_polynomial_eval(pwave, pos);
734  if (!pwave || !isnormal(ltest) || fabs(ltest) < DBL_EPSILON) {
735  char *msg = cpl_sprintf("Wavelength calibration polynomial for slice %d "
736  "of IFU %hhu is not well defined!", islice + 1,
737  ifu);
738  cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT, "%s", msg);
739  cpl_msg_error(__func__, "%s", msg);
740  cpl_free(msg);
741  cpl_polynomial_delete(pwave);
742  cpl_vector_delete(pos);
743  if (getenv("MUSE_EXPERT_USER")) {
744  continue; /* go the next slice */
745  } else { /* give up completely */
746  cpl_table_delete(geopos);
748  return NULL;
749  }
750  } /* if bad wavecal polynomial */
751 
752  /* get the tracing polynomials for this slice */
753  cpl_polynomial **ptrace = muse_trace_table_get_polys_for_slice(aTrace,
754  islice + 1);
755  if (!ptrace) {
756  char *msg = cpl_sprintf("Tracing polynomials for slice %d of IFU %hhu are"
757  " not well defined!", islice + 1, ifu);
758  cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT, "%s", msg);
759  cpl_msg_error(__func__, "%s", msg);
760  cpl_free(msg);
761  cpl_polynomial_delete(pwave);
762  cpl_vector_delete(pos);
763  if (getenv("MUSE_EXPERT_USER")) {
764  continue; /* go the next slice */
765  } else { /* give up completely */
766  cpl_table_delete(geopos);
768  return NULL;
769  }
770  } /* if no tracing polynomials */
771  unsigned offset = muse_pixtable_origin_set_offset(pt,
772  ptrace[MUSE_TRACE_LEFT],
773  ifu, islice + 1);
774 
775  /* within each slice, loop from bottom to top */
776  int j;
777  for (j = 1; j <= ysize; j++) {
778  /* compute slice center for this vertical position and set the *
779  * edge pixels of slice, so that we can loop over those pixels */
780  double x1 = cpl_polynomial_eval_1d(ptrace[MUSE_TRACE_LEFT], j, NULL),
781  x2 = cpl_polynomial_eval_1d(ptrace[MUSE_TRACE_RIGHT], j, NULL),
782  xcenter = (x1 + x2) / 2.,
783  width = x2 - x1;
784  if (!isnormal(x1) || !isnormal(x2) || x1 < 1 || x2 > xsize || x1 > x2) {
785  cpl_msg_warning(__func__, "slice %2d of IFU %hhu: faulty polynomial "
786  "detected at y=%d (borders: %f ... %f)", islice + 1,
787  ifu, j, x1, x2);
788  break; /* skip the rest of this slice */
789  }
790  /* include pixels that have more than half of them inside the slice */
791  int ileft = ceil(x1),
792  iright = floor(x2);
793  cpl_vector_set(pos, 1, j); /* vertical pos. for wavelength evaluation */
794 
795  /* now loop over all pixels of this slice horizontally */
796  int i;
797  for (i = ileft; i <= iright; i++) {
798  cpl_vector_set(pos, 0, i); /* horiz. pos. for wavelength evaluation */
799 
800  /* do the wavelength evaluation early on, to more quickly *
801  * recover from errors */
802  double lambda = cpl_polynomial_eval(pwave, pos);
803  if (lambda < 3000. || lambda > 11000. || !isfinite(lambda)) {
804  continue; /* skip this pixel */
805  }
806 
807  /* Relative horizontal pixel position within slice; has to be scaled *
808  * by the width of the slice on the CCD. Slice width is measured *
809  * horizontally, not along iso-lambda lines, the difference for *
810  * straight lines is below 0.015 pix (0.02%), and only slightly *
811  * larger for curved lines, so this simplification should not have a *
812  * negative effect on anything. We compute (xcenter - i) to account *
813  * for the reversed sky orientation within each slice. */
814  float dx = (xcenter - i) / width
815  * (slice_width ? slice_width[islice] : kMuseSliceNominalWidth);
816 #if DEBUG_PIXTABLE_CREATION /* check how many pixels extend over the nominal width */
817  if (fabs(dx) > 37.5) {
818  cpl_msg_debug(__func__, "%d - %f -> %f", i, xcenter, dx);
819  }
820 #endif
821 
822 #if CREATE_MINIMAL_PIXTABLE
823  /* if we need to create a very small pixel table for debugging */
824  if (lambda < 6495. || lambda > 6505. || fabs(dx) > 10.) {
825  continue;
826  }
827 #endif
828 
829  /* this pixel is set up to be inside the slice, write it to the table */
830  cpl_size irow = itablerow++; /* started at 0, so increment after assignment */
831 #if !PIXTABLE_CREATE_CCDSIZED /* overflow check, and possibly enlarge the table */
832  if (irow + 1 > (cpl_size)muse_pixtable_get_nrow(pt)) {
833 #if DEBUG_PIXTABLE_CREATION /* make sure that we are filling the right part */
834  printf("pixel table before enlargement:\n");
835  cpl_table_dump(pt->table, irow-5, 10, stdout);
836  fflush(stdout);
837 #endif
838  /* If the table is too small, we can as well add a lot more rows *
839  * right now, than adding them in small portions. If we add too many *
840  * this will only marginally decrease the speed of this procedure */
841  const cpl_size nr = 1000000;
842  cpl_msg_debug(__func__, "expand table to %"CPL_SIZE_FORMAT" rows "
843  "(add %"CPL_SIZE_FORMAT" lines)!", irow + nr, nr);
844  cpl_table_set_size(pt->table, irow + nr);
845 
846  /* again fill the new section of each column */
847  cpl_table_fill_column_window_float(pt->table, MUSE_PIXTABLE_XPOS, irow, nr, 0);
848  cpl_table_fill_column_window_float(pt->table, MUSE_PIXTABLE_YPOS, irow, nr, 0);
849  cpl_table_fill_column_window_float(pt->table, MUSE_PIXTABLE_LAMBDA, irow, nr, 0);
850  cpl_table_fill_column_window_float(pt->table, MUSE_PIXTABLE_DATA, irow, nr, 0);
851  cpl_table_fill_column_window_float(pt->table, MUSE_PIXTABLE_STAT, irow, nr, 0);
852  cpl_table_fill_column_window_int(pt->table, MUSE_PIXTABLE_DQ, irow, nr, 0);
853  cpl_table_fill_column_window_int(pt->table, MUSE_PIXTABLE_ORIGIN, irow, nr, 0);
854 
855  /* re-get columns as pointers, the memory location might have changed */
856  cdata_xpos = cpl_table_get_data_float(pt->table, MUSE_PIXTABLE_XPOS);
857  cdata_ypos = cpl_table_get_data_float(pt->table, MUSE_PIXTABLE_YPOS);
858  cdata_lambda = cpl_table_get_data_float(pt->table, MUSE_PIXTABLE_LAMBDA);
859  cdata_data = cpl_table_get_data_float(pt->table, MUSE_PIXTABLE_DATA);
860  cdata_stat = cpl_table_get_data_float(pt->table, MUSE_PIXTABLE_STAT);
861  cdata_dq = cpl_table_get_data_int(pt->table, MUSE_PIXTABLE_DQ);
862  cdata_origin = (uint32_t *)cpl_table_get_data_int(pt->table,
864 #if DEBUG_PIXTABLE_CREATION /* make sure that we have filled the right part */
865  printf("pixel table after filling new rows:\n");
866  cpl_table_dump(pt->table, irow-5, 10, stdout);
867  fflush(stdout);
868 #endif
869  } /* if new rows needed */
870 #endif /* !PIXTABLE_CREATE_CCDSIZED */
871 
872  /* Assign the coordinates (slice angle is in degrees). */
873  double angle = slice_angle ? slice_angle[islice] * CPL_MATH_RAD_DEG
874  : 0.;
875  cdata_xpos[irow] = aGeoTable
876  ? dx * cos(angle) + slice_x[islice]
877  : dx + islice * kMuseSliceHiLikelyWidth
878  + kMuseSliceHiLikelyWidth / 2.;
879  cdata_ypos[irow] = aGeoTable
880  ? dx * sin(angle) + slice_y[islice]
881  : 0.; /* no vertical position for simple case */
882 #if DEBUG_PIXTABLE_CREATION
883  if (!isfinite(cdata_xpos[irow]) ||
884  cdata_xpos[irow] < -200 || cdata_xpos[irow] > 4000) {
885  cpl_msg_debug(__func__, "weird data in x: %e %e",
886  cdata_xpos[irow], cdata_ypos[irow]);
887  }
888  if (!isfinite(cdata_ypos[irow]) ||
889  cdata_ypos[irow] < -200 || cdata_ypos[irow] > 200) {
890  cpl_msg_debug(__func__, "weird data in y: %e %e",
891  cdata_xpos[irow], cdata_ypos[irow]);
892  }
893 #endif
894  cdata_lambda[irow] = lambda;
895 
896  /* assign the data values */
897  cdata_data[irow] = pixdata[(i-1) + (j-1)*xsize];
898  cdata_dq[irow] = pixdq[(i-1) + (j-1)*xsize];
899  cdata_stat[irow] = pixstat[(i-1) + (j-1)*xsize];
900  /* assign origin and slice number */
901  cdata_origin[irow] = muse_pixtable_origin_encode_fast(i, j, ifu,
902  islice + 1,
903  offset);
904  } /* for i (horizontal pixels) */
905  } /* for j (vertical direction) */
906 
907  /* we are now done with this slice, clean up */
908  muse_trace_polys_delete(ptrace);
909  cpl_polynomial_delete(pwave);
910  cpl_vector_delete(pos);
911  } /* for islice */
912  cpl_table_delete(geopos);
913  cpl_msg_debug(__func__, "IFU %hhu took %gs (CPU time), %gs (wall-clock)",
914  ifu, cpl_test_get_cputime() - cputime,
915  cpl_test_get_walltime() - walltime);
916 
917  /* itablerow is the last written to location (starting at zero), so add one */
918  if (muse_pixtable_get_nrow(pt) > itablerow) {
919  cpl_msg_debug(__func__, "Trimming pixel table of IFU %hhu to %"CPL_SIZE_FORMAT
920  " of %"CPL_SIZE_FORMAT" rows", ifu, itablerow,
922 #if DEBUG_PIXTABLE_CREATION /* make sure that we are trimming correctly */
923  printf("end of used part of pixel table before trimming:\n");
924  cpl_table_dump(pt->table, itablerow-5, 10, stdout);
925  fflush(stdout);
926 #endif
927  cpl_table_set_size(pt->table, itablerow);
928  }
929  /* now that we are done creating it, we can also *
930  * compute its extremes for the first time */
932 #if DEBUG_PIXTABLE_CREATION
933  printf("beginning of pixel table before returning:\n");
934  cpl_table_dump(pt->table, 0, 5, stdout);
935  fflush(stdout);
936  printf("end of pixel table before returning:\n");
937  cpl_table_dump(pt->table, muse_pixtable_get_nrow(pt)-5, 10, stdout);
938  fflush(stdout);
939 #endif
940 
941  return pt;
942 } /* muse_pixtable_create() */
943 
944 /*---------------------------------------------------------------------------*/
957 /*---------------------------------------------------------------------------*/
960 {
961  if (aPixtable == NULL) {
962  return NULL;
963  }
964  muse_pixtable *pt = cpl_calloc(1, sizeof(muse_pixtable));
965  pt->table = cpl_table_duplicate(aPixtable->table);
966  pt->header = cpl_propertylist_duplicate(aPixtable->header);
967  return pt;
968 }
969 
970 /*---------------------------------------------------------------------------*/
979 /*---------------------------------------------------------------------------*/
980 void
982 {
983  /* if we already get NULL, we don't need to do anything */
984  if (!aPixtable) {
985  return;
986  }
987 
988  /* delete the table portion */
989  cpl_table_delete(aPixtable->table);
990  aPixtable->table = NULL;
991 
992  /* delete the header portion */
993  cpl_propertylist_delete(aPixtable->header);
994  aPixtable->header = NULL;
995 
996  /* delete the flat-field spectrum */
997  cpl_table_delete(aPixtable->ffspec);
998  aPixtable->ffspec = NULL;
999 
1000  cpl_free(aPixtable);
1001 } /* muse_pixtable_delete() */
1002 
1003 /*---------------------------------------------------------------------------*/
1028 /*---------------------------------------------------------------------------*/
1029 cpl_error_code
1031  cpl_table *aTrace, cpl_table *aWave, float aSampling)
1032 {
1033  cpl_ensure_code(aPixtable, CPL_ERROR_NULL_INPUT);
1034 
1035  /* discard a pre-existing flat-field spectrum */
1036  if (aPixtable->ffspec) {
1037  cpl_table_delete(aPixtable->ffspec);
1038  }
1039 
1040  /* create pixel table, without geometry table is fine */
1041  muse_pixtable *ptflat = muse_pixtable_create(aFF, aTrace, aWave, NULL);
1042  if (!ptflat) {
1043  return cpl_error_get_code();
1044  }
1045  aPixtable->ffspec = muse_resampling_spectrum(ptflat, aSampling);
1046  muse_pixtable_delete(ptflat);
1047  /* there does not seem to be a way to generate a muse_resampling_spectrum() *
1048  * failure other than with a NULL pixel table, so do not test for that */
1049 
1050  /* muse_resampling_spectrum() also creates variance and bad pixel *
1051  * spectra in extra columns which we do not need to carry around here */
1052  cpl_table_erase_column(aPixtable->ffspec, "stat");
1053  cpl_table_erase_column(aPixtable->ffspec, "dq");
1054 
1055  return CPL_ERROR_NONE;
1056 } /* muse_pixtable_append_ff() */
1057 
1058 /*---------------------------------------------------------------------------*/
1070 /*---------------------------------------------------------------------------*/
1071 static cpl_error_code
1072 muse_pixtable_save_ffspec(muse_pixtable *aPixtable, const char *aFilename)
1073 {
1074  cpl_ensure_code(aPixtable, CPL_ERROR_NULL_INPUT);
1075 
1076  /* if the flat-field spectrum does not exist, return early */
1077  if (!aPixtable->ffspec) {
1078  /* but warn, if the flag was set in the primary header *
1079  * (too late now to do anything about it) */
1080  return CPL_ERROR_NONE;
1081  }
1082 
1083  /* save the flat-field in a new file extension */
1084  cpl_propertylist *hext = cpl_propertylist_new();
1085  cpl_propertylist_append_string(hext, "EXTNAME", MUSE_PIXTABLE_FF_EXT);
1086  cpl_error_code rc = cpl_table_save(aPixtable->ffspec, NULL, hext, aFilename,
1087  CPL_IO_EXTEND);
1088  cpl_propertylist_delete(hext);
1089 
1090  return rc;
1091 } /* muse_pixtable_save_ffspec() */
1092 
1093 /*---------------------------------------------------------------------------*/
1114 /*---------------------------------------------------------------------------*/
1115 static cpl_error_code
1116 muse_pixtable_save_image(muse_pixtable *aPixtable, const char *aFilename)
1117 {
1118  const char *id = "muse_pixtable_save"; /* pretend to be in that function */
1119  cpl_ensure_code(aPixtable, CPL_ERROR_NULL_INPUT);
1120  cpl_size nrow = muse_pixtable_get_nrow(aPixtable);
1121  cpl_ensure_code(nrow > 0, CPL_ERROR_ILLEGAL_INPUT);
1122 
1123  cpl_errorstate state = cpl_errorstate_get();
1124  cpl_array *columns = cpl_table_get_column_names(aPixtable->table);
1125  int icol, ncol = cpl_array_get_size(columns);
1126  for (icol = 0; icol < ncol; icol++) {
1127  const char *colname = cpl_array_get_string(columns, icol);
1128  /* access data buffer and wrap as corresponding CPL image type */
1129  cpl_type type = cpl_table_get_column_type(aPixtable->table, colname);
1130  cpl_image *image = NULL;
1131  switch (type) {
1132  case CPL_TYPE_FLOAT: {
1133  float *data = cpl_table_get_data_float(aPixtable->table, colname);
1134  image = cpl_image_wrap_float(1, nrow, data);
1135  break;
1136  }
1137  case CPL_TYPE_INT: {
1138  int *data = cpl_table_get_data_int(aPixtable->table, colname);
1139  image = cpl_image_wrap_int(1, nrow, data);
1140  break;
1141  }
1142  default:
1143  cpl_error_set_message(id, CPL_ERROR_UNSUPPORTED_MODE, "type \"%s\" (of "
1144  "column %s) is not supported for MUSE pixel tables",
1145  cpl_type_get_name(type), colname);
1146  continue;
1147  } /* switch */
1148  /* now construct minimal extension header with EXTNAME and *
1149  * BUNIT and save the image with its internal data type */
1150  cpl_propertylist *hext = cpl_propertylist_new();
1151  cpl_propertylist_append_string(hext, "EXTNAME", colname);
1152  const char *unit = cpl_table_get_column_unit(aPixtable->table, colname);
1153  if (unit) {
1154  cpl_propertylist_append_string(hext, "BUNIT", unit);
1155  }
1156  cpl_image_save(image, aFilename, CPL_TYPE_UNSPECIFIED, hext, CPL_IO_EXTEND);
1157  cpl_image_unwrap(image);
1158  cpl_propertylist_delete(hext);
1159  } /* for icol (all table columns) */
1160  cpl_array_delete(columns);
1161 
1162  /* append the flat-field spectrum to the file (if it exists) */
1163  muse_pixtable_save_ffspec(aPixtable, aFilename);
1164 
1165  return cpl_errorstate_is_equal(state) ? CPL_ERROR_NONE : cpl_error_get_code();
1166 } /* muse_pixtable_save_image() */
1167 
1168 /*---------------------------------------------------------------------------*/
1189 /*---------------------------------------------------------------------------*/
1190 cpl_error_code
1191 muse_pixtable_save(muse_pixtable *aPixtable, const char *aFilename)
1192 {
1193  cpl_ensure_code(aPixtable, CPL_ERROR_NULL_INPUT);
1194 
1195  /* let CPL do the real work and propagate any errors and errorstates that *
1196  * it generates; the secondary header should be generated automatically *
1197  * from the table data we have set in muse_pixtable_create(). */
1198  /* save headers separately to keep WCS information */
1199  cpl_error_code rc = cpl_propertylist_save(aPixtable->header, aFilename,
1200  CPL_IO_CREATE);
1201  if (rc != CPL_ERROR_NONE) {
1202  cpl_error_set_message(__func__, rc, "could not save FITS header of pixel "
1203  "table \"%s\"", aFilename);
1204  return rc;
1205  }
1206  /* test if we need to save it as image */
1207  cpl_boolean astable = getenv("MUSE_PIXTABLE_SAVE_AS_TABLE")
1208  && atoi(getenv("MUSE_PIXTABLE_SAVE_AS_TABLE")) > 0;
1209  if (astable) {
1210  cpl_msg_debug(__func__, "Saving pixel table \"%s\" as binary table",
1211  aFilename);
1212  rc = cpl_table_save(aPixtable->table, NULL, NULL, aFilename, CPL_IO_EXTEND);
1213  /* append the flat-field spectrum to the file (if it exists) */
1214  muse_pixtable_save_ffspec(aPixtable, aFilename);
1215  } else {
1216  rc = muse_pixtable_save_image(aPixtable, aFilename);
1217  }
1218  return rc;
1219 } /* muse_pixtable_save() */
1220 
1221 /*---------------------------------------------------------------------------*/
1246 /*---------------------------------------------------------------------------*/
1247 static cpl_table *
1248 muse_pixtable_load_window_image(const char *aFilename,
1249  cpl_size aStart, cpl_size aNRows)
1250 {
1251  const char *id = "muse_pixtable_load"; /* pretend to be in that function */
1252  /* determine "image" window to read, the y-direction *
1253  * in the image is the dimension of the rows */
1254  cpl_size y1 = aStart + 1,
1255  y2 = aStart + aNRows;
1256  /* get header of "data" extension, and check/correct the size of the dataset */
1257  int ext = cpl_fits_find_extension(aFilename, MUSE_PIXTABLE_DATA);
1258  cpl_propertylist *hdata = cpl_propertylist_load(aFilename, ext);
1259  /* get NAXIS2 which contains contains the number of rows */
1260  cpl_size naxis2 = muse_pfits_get_naxis(hdata, 2);
1261  if (y2 > naxis2) {
1262  y2 = naxis2;
1263  }
1264  cpl_propertylist_delete(hdata);
1265 
1266  /* find and load all image extensions using their EXTNAMEs */
1267  cpl_size nrow = 0;
1268  cpl_table *table = cpl_table_new(nrow);
1269  int iext, next = cpl_fits_count_extensions(aFilename);
1270  for (iext = 1; iext <= next; iext++) {
1271  /* first check, that we are not loading the flat-field extension */
1272  cpl_propertylist *hext = cpl_propertylist_load(aFilename, iext);
1273  const char *colname = muse_pfits_get_extname(hext);
1274  if (!strncmp(colname, MUSE_PIXTABLE_FF_EXT, strlen(MUSE_PIXTABLE_FF_EXT) + 1)) {
1275  /* this will be loaded later */
1276  cpl_propertylist_delete(hext);
1277  continue;
1278  }
1279 
1280  /* now load the column from the image extension */
1281  cpl_errorstate ps = cpl_errorstate_get();
1282  cpl_image *column = cpl_image_load_window(aFilename, CPL_TYPE_UNSPECIFIED,
1283  0, iext, 1, y1, 1, y2);
1284  if (!column || !cpl_errorstate_is_equal(ps)) {
1285  cpl_image_delete(column);
1286  cpl_error_set_message(id, cpl_error_get_code(), "could not load extension"
1287  " %d of pixel table \"%s\"", iext, aFilename);
1288  cpl_propertylist_delete(hext);
1289  continue;
1290  }
1291  cpl_size nrows = cpl_image_get_size_x(column)
1292  * cpl_image_get_size_y(column);
1293  if (nrow < 1) {
1294  cpl_table_set_size(table, nrows);
1295  nrow = nrows;
1296  } else if (nrows != nrow) {
1297  cpl_error_set_message(id, CPL_ERROR_INCOMPATIBLE_INPUT, "size of column "
1298  "%s does not match", colname);
1299  cpl_propertylist_delete(hext);
1300  cpl_image_delete(column);
1301  continue;
1302  }
1303  cpl_type type = cpl_image_get_type(column);
1304  switch (type) {
1305  case CPL_TYPE_FLOAT:
1306  cpl_table_wrap_float(table, cpl_image_unwrap(column), colname);
1307  break;
1308  case CPL_TYPE_INT:
1309  cpl_table_wrap_int(table, cpl_image_unwrap(column), colname);
1310  break;
1311  default:
1312  cpl_error_set_message(id, CPL_ERROR_UNSUPPORTED_MODE, "type \"%s\" (of "
1313  "column %s) is not supported for MUSE pixel tables",
1314  cpl_type_get_name(type), colname);
1315  } /* switch */
1316  cpl_errorstate state = cpl_errorstate_get();
1317  const char *unit = muse_pfits_get_bunit(hext);
1318  if (!cpl_errorstate_is_equal(state)) {
1319  cpl_errorstate_set(state);
1320  }
1321  if (unit) {
1322  cpl_table_set_column_unit(table, colname, unit);
1323  }
1324  cpl_propertylist_delete(hext);
1325  } /* for iext (all FITS extensions) */
1326  return table;
1327 } /* muse_pixtable_load_window_image() */
1328 
1329 /*---------------------------------------------------------------------------*/
1344 /*---------------------------------------------------------------------------*/
1345 static cpl_error_code
1346 muse_pixtable_load_ffspec(const char *aFilename, muse_pixtable *aPixtable)
1347 {
1348  const char *id = "muse_pixtable_load"; /* pretend to be in that function */
1349  cpl_ensure_code(aFilename && aPixtable, CPL_ERROR_NULL_INPUT);
1350 
1351  int iext = cpl_fits_find_extension(aFilename, MUSE_PIXTABLE_FF_EXT);
1352  /* if the flat-field extension was not found, just quietly *
1353  * return -- or make some noise if the flag was set */
1354  if (iext <= 0) {
1355  return CPL_ERROR_NONE; /* quietly return */
1356  }
1357 
1358  cpl_errorstate es = cpl_errorstate_get();
1359  aPixtable->ffspec = cpl_table_load(aFilename, iext, 1);
1360  if (!cpl_errorstate_is_equal(es)) {
1361  cpl_msg_warning(id, "Pixel table flat-field spectrum extension %s "
1362  "exists in \"%s\", but cannot be loaded!",
1363  MUSE_PIXTABLE_FF_EXT, aFilename);
1364 
1365  /* clean up again, in case something (broken?) was loaded */
1366  cpl_table_delete(aPixtable->ffspec);
1367  aPixtable->ffspec = NULL;
1368  /* swallow the error code, since the missing flat-field only becomes *
1369  * critical, if this pixel table is later combined with others */
1370  cpl_errorstate_set(es);
1371  }
1372 
1373  return CPL_ERROR_NONE;
1374 } /* muse_pixtable_load_ffspec() */
1375 
1376 /*---------------------------------------------------------------------------*/
1401 /*---------------------------------------------------------------------------*/
1402 muse_pixtable *
1403 muse_pixtable_load_window(const char *aFilename,
1404  cpl_size aStart, cpl_size aNRows)
1405 {
1406  /* create the output pixel table structure, and NULL out the components */
1407  muse_pixtable *pt = cpl_calloc(1, sizeof(muse_pixtable));
1408 
1409  /* load the header and return with an error if something went wrong */
1410  cpl_errorstate prestate = cpl_errorstate_get();
1411  pt->header = cpl_propertylist_load(aFilename, 0);
1412  cpl_ensure(cpl_errorstate_is_equal(prestate) && pt->header,
1413  cpl_error_get_code(), NULL);
1415  cpl_msg_error(__func__, "unknown pixel table type found in \"%s\"", aFilename);
1417  return NULL;
1418  }
1419 
1420  /* determine the type of data (FITS binary table extension or *
1421  * FITS multi-extension images) by looking at the first extension */
1422  cpl_propertylist *hext = cpl_propertylist_load(aFilename, 1);
1423  cpl_boolean asimage = !strcmp(cpl_propertylist_get_string(hext, "XTENSION"),
1424  "IMAGE");
1425  cpl_propertylist_delete(hext);
1426 
1427  /* load the table itself and return an error if something went wrong */
1428  if (asimage) {
1429  cpl_msg_info(__func__, "Loading pixel table \"%s\" (image format)",
1430  aFilename);
1431  pt->table = muse_pixtable_load_window_image(aFilename, aStart, aNRows);
1432  } else {
1433  cpl_msg_info(__func__, "Loading pixel table \"%s\" (bintable format)",
1434  aFilename);
1435  pt->table = cpl_table_load_window(aFilename, 1, 0, NULL, aStart, aNRows);
1436  }
1437  if (!cpl_errorstate_is_equal(prestate) || !pt->table) {
1438  cpl_msg_error(__func__, "Failed to load table part of pixel table \"%s\"",
1439  aFilename);
1441  return NULL;
1442  }
1443  cpl_error_code rc = muse_cpltable_check(pt->table, muse_pixtable_def);
1444  if (rc != CPL_ERROR_NONE) {
1445  cpl_error_set_message(__func__, rc, "pixel table \"%s\" does not contain "
1446  "all expected columns", aFilename);
1447  }
1448 
1449  /* load the flat-field spectrum, if it exists; same for both cases */
1450  muse_pixtable_load_ffspec(aFilename, pt);
1451  return pt;
1452 } /* muse_pixtable_load_window() */
1453 
1454 /*---------------------------------------------------------------------------*/
1474 /*---------------------------------------------------------------------------*/
1475 muse_pixtable *
1476 muse_pixtable_load(const char *aFilename)
1477 {
1478  /* find out original table length from the header of the FITS binary table */
1479  cpl_errorstate prestate = cpl_errorstate_get();
1480  cpl_propertylist *theader = cpl_propertylist_load(aFilename, 1);
1481  cpl_ensure(cpl_errorstate_is_equal(prestate) && theader,
1482  cpl_error_get_code(), NULL);
1483  cpl_size nrow = cpl_propertylist_get_long_long(theader, "NAXIS2");
1484  cpl_propertylist_delete(theader);
1485 
1486  /* now that we know the full size, load all rows */
1487  muse_pixtable *pt = muse_pixtable_load_window(aFilename, 0, nrow);
1488  return pt;
1489 } /* muse_pixtable_load() */
1490 
1491 /*---------------------------------------------------------------------------*/
1508 /*---------------------------------------------------------------------------*/
1509 muse_pixtable *
1511  double aLambdaMin, double aLambdaMax)
1512 {
1513  muse_pixtable *pt = muse_pixtable_load(aFilename);
1514  if (!pt) {
1515  return NULL;
1516  }
1517  cpl_error_code rc = muse_pixtable_restrict_wavelength(pt, aLambdaMin,
1518  aLambdaMax);
1519  if (rc != CPL_ERROR_NONE) {
1521  return NULL;
1522  }
1523  if (muse_pixtable_get_nrow(pt) < 1) {
1524  cpl_msg_error(__func__, "Pixel table contains no entries after cutting to "
1525  "%.3f..%.3f Angstrom", aLambdaMin, aLambdaMax);
1526  cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
1528  return NULL;
1529  }
1530  return pt;
1531 } /* muse_pixtable_load_restricted_wavelength() */
1532 
1533 /*---------------------------------------------------------------------------*/
1580 /*---------------------------------------------------------------------------*/
1581 muse_pixtable *
1582 muse_pixtable_load_merge_channels(cpl_table *aExposureList,
1583  double aLambdaMin, double aLambdaMax)
1584 {
1585  cpl_ensure(aExposureList, CPL_ERROR_NULL_INPUT, NULL);
1586 
1587  muse_pixtable *pt = NULL; /* the big output pixel table */
1588  if (cpl_table_has_column(aExposureList, "00")) {
1589  const char *filename = cpl_table_get_string(aExposureList, "00", 0);
1590  if (filename) {
1591  pt = muse_pixtable_load_restricted_wavelength(filename, aLambdaMin,
1592  aLambdaMax);
1593  }
1594  if (pt) {
1595  return pt;
1596  }
1597  } /* if table with "00" column */
1598 
1599  /* create array to hold interpolation points for flat-field spectrum */
1600  cpl_array *lambdas = cpl_array_new((kMuseLambdaMaxX - kMuseLambdaMinX)
1601  / kMuseSpectralSamplingA + 1,
1602  CPL_TYPE_DOUBLE),
1603  *ffspec = NULL; /* for output average flat-field spectrum data */
1604  int j, nwl = cpl_array_get_size(lambdas);
1605  for (j = 0; j < nwl; j++) {
1606  cpl_array_set_double(lambdas, j, kMuseLambdaMinX
1607  + j * kMuseSpectralSamplingA);
1608  }
1609 
1610  cpl_boolean isfirst = CPL_TRUE;
1611  double fluxlref = 0, /* propagated from lamp flat */
1612  fluxsref = 0; /* propagated from sky flat */
1613  int i, nifu = 0, nflats = 0;
1614  for (i = 1; i <= kMuseNumIFUs; i++) {
1615  char *colname = cpl_sprintf("%02d", i);
1616  const char *filename = cpl_table_get_string(aExposureList, colname, 0);
1617  cpl_free(colname);
1618  if (!filename) {
1619  cpl_msg_warning(__func__, "Channel for IFU %02d is missing", i);
1620  continue;
1621  }
1623  aLambdaMin,
1624  aLambdaMax);
1625  if (!onept) {
1626  cpl_msg_error(__func__, "failed to load pixel table from \"%s\"",
1627  filename);
1628  cpl_array_delete(lambdas);
1629  return pt;
1630  }
1631  nifu++;
1632 
1633  /* if this is the first pixel table, use its data to create the big one */
1634  if (isfirst) {
1635  pt = onept;
1636  cpl_msg_debug(__func__, "loaded pixel table with %"CPL_SIZE_FORMAT" rows",
1638  isfirst = CPL_FALSE;
1639  cpl_errorstate prestate = cpl_errorstate_get();
1640  fluxsref = cpl_propertylist_get_double(pt->header, MUSE_HDR_FLAT_FLUX_SKY);
1641  fluxlref = cpl_propertylist_get_double(pt->header, MUSE_HDR_FLAT_FLUX_LAMP);
1642  if (fluxsref == 0. && fluxlref == 0. && !cpl_errorstate_is_equal(prestate)) {
1643  /* If the flat headers are both missing, this can only mean that *
1644  * the exposure was previously merged into the pixel table we just *
1645  * loaded. Then we can recover from the error and stop right here. */
1646  cpl_msg_debug(__func__, "\"%s\" was previously merged (got \"%s\" when"
1647  " asking for flat-field fluxes)", filename,
1648  cpl_error_get_message());
1649  cpl_errorstate_set(prestate);
1650  break;
1651  }
1652  if (fluxsref == 0. && fluxlref > 0. && !cpl_errorstate_is_equal(prestate)) {
1653  /* only the sky-flat level was missing, output a warning */
1654  cpl_msg_warning(__func__, "only found reference lamp-flat flux (%e) in "
1655  "\"%s\", flux levels may vary between IFUs!", fluxlref,
1656  filename);
1657  cpl_errorstate_set(prestate);
1658  } else {
1659  cpl_msg_debug(__func__, "reference flat fluxes sky: %e lamp: %e",
1660  fluxsref, fluxlref);
1661  }
1662  cpl_propertylist_erase(pt->header, MUSE_HDR_FLAT_FLUX_SKY);
1663  cpl_propertylist_erase(pt->header, MUSE_HDR_FLAT_FLUX_LAMP);
1664 
1665  /* for the first table, directly resample the flat-field spectrum */
1666  if (pt->ffspec) { /* check, to keep backward compatibility */
1667  ffspec = muse_cplarray_interpolate_table_linear(lambdas, pt->ffspec,
1668  "lambda", "data");
1669  /* delete the flat-field spectrum component, since *
1670  * this will be the flat-field corrected output table */
1671  cpl_table_delete(pt->ffspec);
1672  pt->ffspec = NULL;
1673  nflats++;
1674  }
1675  continue;
1676  } /* if isfirst */
1677 
1678  /* copy the offset headers to the output pixel table */
1679  muse_pixtable_origin_copy_offsets(pt, onept, 0);
1680 
1681  /* to fix the flux offset (using propagated flat-field fluxes) */
1682  cpl_errorstate state = cpl_errorstate_get();
1683  double fluxs = cpl_propertylist_get_double(onept->header,
1684  MUSE_HDR_FLAT_FLUX_SKY),
1685  fluxl = cpl_propertylist_get_double(onept->header,
1686  MUSE_HDR_FLAT_FLUX_LAMP),
1687  scale = 1.;
1688  if (fluxsref > 0. && fluxs > 0.) { /* prefer sky flat flux scaling */
1689  scale = fluxs / fluxsref;
1690  } else if (fluxlref > 0. && fluxl > 0.) {
1691  scale = fluxl / fluxlref;
1692  if (!cpl_errorstate_is_equal(state)) {
1693  cpl_msg_warning(__func__, "only found relative lamp-flat flux (%e) in "
1694  "\"%s\", flux levels may vary between IFUs!", fluxl,
1695  filename);
1696  cpl_errorstate_set(state);
1697  } /* if bad error state */
1698  } /* else: use lamp-flat values */
1699  muse_pixtable_flux_multiply(onept, 1. / scale); /* inverse of scale here */
1700 
1701  /* resample ff-spectrum onto the regular sampling, given *
1702  * by the lambda vector created at the beginning */
1703  if (onept->ffspec) { /* check, to keep backward compatibility */
1704  cpl_array *ffout = muse_cplarray_interpolate_table_linear(lambdas,
1705  onept->ffspec,
1706  "lambda", "data");
1707  if (ffspec) { /* if the first one had a flat-field spectrum */
1708  cpl_array_add(ffspec, ffout); /* add up to average the resulting spectrum */
1709  }
1710  cpl_array_delete(ffout);
1711  nflats++;
1712  } /* if ff-spectrum */
1713 
1714  /* append this pixel table to the big one */
1715  cpl_table_insert(pt->table, onept->table, muse_pixtable_get_nrow(pt));
1716  cpl_msg_debug(__func__, "big pixel table now has %"CPL_SIZE_FORMAT" entries,"
1717  " scale was %e (flat fluxes sky: %e lamp: %e)",
1718  muse_pixtable_get_nrow(pt), scale, fluxs, fluxl);
1719 
1720  /* the insertion duplicates the data, so we can delete it now */
1721  muse_pixtable_delete(onept);
1722  } /* for i (all IFUs) */
1723 
1724  /* If some but not all pixel tables came with the flat-field spectrum, *
1725  * we may get unforeseen results when doing flux calibration etc., so *
1726  * better fail completely! */
1727  if (nflats > 0 && nifu != nflats) {
1728  cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT, "Only %d of %d "
1729  "pixel tables of this exposure came with a flat-field"
1730  " spectrum, cannot continue!", nflats, nifu);
1731  cpl_array_delete(lambdas);
1732  cpl_array_delete(ffspec);
1734  return NULL;
1735  } /* if nflats */
1736 
1737  if (ffspec) {
1738  /* we need the average flat-field spectrum, so *
1739  * divide by the number of IFUs that were combined */
1740  cpl_array_divide_scalar(ffspec, nflats);
1741  cpl_msg_debug(__func__, "Average of flat-field spectrum: %.4f",
1742  cpl_array_get_mean(ffspec));
1743 
1744  /* now correct the output pixel table by the averaged flat-field spectrum */
1745  muse_pixtable_spectrum_apply(pt, lambdas, ffspec,
1746  MUSE_PIXTABLE_OPERATION_MULTIPLY);
1747 
1748  /* mark the header as having gotten the correction */
1749  cpl_propertylist_update_int(pt->header, MUSE_HDR_PT_FFCORR, nflats);
1750  cpl_propertylist_set_comment(pt->header, MUSE_HDR_PT_FFCORR,
1751  MUSE_HDR_PT_FFCORR_COMMENT);
1752 
1753  /* still keep the (rather small) flat-field spectrum around, *
1754  * at least for debugging */
1755  pt->ffspec = cpl_table_new(cpl_array_get_size(lambdas));
1756  cpl_table_new_column(pt->ffspec, MUSE_PIXTABLE_FFLAMBDA, CPL_TYPE_DOUBLE);
1757  cpl_table_new_column(pt->ffspec, MUSE_PIXTABLE_FFDATA, CPL_TYPE_DOUBLE);
1760  cpl_array_delete(ffspec);
1761  /* everything before was done on a fixed grid to encompass *
1762  * all possible MUSE data, now we have to get rid of the invalid *
1763  * entries at the extreme ends of the flat-field spectrum */
1764  cpl_table_erase_invalid(pt->ffspec);
1765  } /* if ffspec */
1766  cpl_array_delete(lambdas);
1767 
1768  /* now that we merged all pixel tables, compute the full extent */
1770  if (!pt) {
1771  cpl_error_set_message(__func__, CPL_ERROR_FILE_NOT_FOUND,
1772  "None of the pixel tables could be loaded");
1773  return NULL; /* exit here to not crash below! */
1774  }
1775  /* data does not belong to any one IFU any more, so *
1776  * remove the EXTNAME header that contains CHAN%02d */
1777  cpl_propertylist_erase_regexp(pt->header, "^EXTNAME|"MUSE_HDR_PT_ILLUM_REGEXP,
1778  0);
1779  /* also remove chip-specific info that doesn't belong to the merged data */
1780  cpl_propertylist_erase_regexp(pt->header, "ESO DET (CHIP|OUT) ", 0);
1781  cpl_propertylist_erase_regexp(pt->header, "ESO DET2 ", 0);
1782  /* add the status header */
1783  cpl_propertylist_update_int(pt->header, MUSE_HDR_PT_MERGED, nifu);
1784  cpl_propertylist_set_comment(pt->header, MUSE_HDR_PT_MERGED,
1785  MUSE_HDR_PT_MERGED_COMMENT);
1786  return pt;
1787 } /* muse_pixtable_load_merge_channels() */
1788 
1789 /*---------------------------------------------------------------------------*/
1801 /*---------------------------------------------------------------------------*/
1802 int
1804 {
1805  cpl_ensure(aPixtable, CPL_ERROR_NULL_INPUT, MUSE_PIXTABLE_TYPE_UNKNOWN);
1806  const char *type = cpl_propertylist_get_string(aPixtable->header,
1808 #if 0
1809  cpl_msg_debug(__func__, "pixel table type \"%s\"", type);
1810 #endif
1811  if (!type) {
1813  }
1814  if (!strncmp(type, MUSE_PIXTABLE_STRING_FULL,
1815  strlen(MUSE_PIXTABLE_STRING_FULL) + 1)) {
1816  return MUSE_PIXTABLE_TYPE_FULL;
1817  } else if (!strncmp(type, MUSE_PIXTABLE_STRING_SIMPLE,
1818  strlen(MUSE_PIXTABLE_STRING_SIMPLE) + 1)) {
1820  }
1821 
1823 } /* muse_pixtable_get_type() */
1824 
1825 /*---------------------------------------------------------------------------*/
1834 /*---------------------------------------------------------------------------*/
1835 cpl_size
1837 {
1838  /* make sure that we are getting valid input */
1839  cpl_ensure(aPixtable, CPL_ERROR_NULL_INPUT, 0);
1840  cpl_ensure(aPixtable->table, CPL_ERROR_NULL_INPUT, 0);
1841 
1842  /* let CPL do the work */
1843  return cpl_table_get_nrow(aPixtable->table);
1844 } /* muse_pixtable_get_nrow() */
1845 
1846 /*---------------------------------------------------------------------------*/
1862 /*---------------------------------------------------------------------------*/
1863 cpl_error_code
1865 {
1866  /* make sure that we are getting valid input */
1867  cpl_ensure_code(aPixtable && aPixtable->table && aPixtable->header,
1868  CPL_ERROR_NULL_INPUT);
1869  cpl_ensure_code(muse_cpltable_check(aPixtable->table, muse_pixtable_def)
1870  == CPL_ERROR_NONE, CPL_ERROR_DATA_NOT_FOUND);
1871 
1872  if (muse_pixtable_get_nrow(aPixtable) == 0) {
1873  return CPL_ERROR_NONE;
1874  }
1875 
1876  float *cdata_xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
1877  *cdata_ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS),
1878  *cdata_lambda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA);
1879  uint32_t *origin = (uint32_t *)cpl_table_get_data_int(aPixtable->table,
1881  /* set to extreme values for a start */
1882  float xlo = FLT_MAX, xhi = -FLT_MAX,
1883  ylo = FLT_MAX, yhi = -FLT_MAX,
1884  llo = FLT_MAX, lhi = -FLT_MAX;
1885  int ifulo = INT_MAX, ifuhi = 0,
1886  slicelo = INT_MAX, slicehi = 0;
1887  cpl_size i, nrow = muse_pixtable_get_nrow(aPixtable);
1888  for (i = 0; i < nrow; i++) {
1889  if (cdata_xpos[i] > xhi) xhi = cdata_xpos[i];
1890  if (cdata_xpos[i] < xlo) xlo = cdata_xpos[i];
1891  if (cdata_ypos[i] > yhi) yhi = cdata_ypos[i];
1892  if (cdata_ypos[i] < ylo) ylo = cdata_ypos[i];
1893  if (cdata_lambda[i] > lhi) lhi = cdata_lambda[i];
1894  if (cdata_lambda[i] < llo) llo = cdata_lambda[i];
1895  int ifu = muse_pixtable_origin_get_ifu_fast(origin[i]),
1896  slice = muse_pixtable_origin_get_slice_fast(origin[i]);
1897  if (ifu > ifuhi) ifuhi = ifu;
1898  if (ifu < ifulo) ifulo = ifu;
1899  if (slice > slicehi) slicehi = slice;
1900  if (slice < slicelo) slicelo = slice;
1901  } /* for i */
1902  char *dodebug = getenv("MUSE_DEBUG_PIXTABLE_LIMITS");
1903  if (dodebug && atoi(dodebug)) {
1904  cpl_msg_debug(__func__, "x: %f...%f, y: %f...%f, lambda: %f...%f, "
1905  "ifu: %d...%d, slice: %d...%d",
1906  xlo, xhi, ylo, yhi, llo, lhi, ifulo, ifuhi, slicelo, slicehi);
1907  }
1908  cpl_propertylist_erase_regexp(aPixtable->header, MUSE_HDR_PT_LIMITS_REGEXP, 0);
1909  double ptxoff = 0., /* zero by default ... */
1910  ptyoff = 0.; /* for pixel coordinates */
1912  ptxoff = muse_pfits_get_crval(aPixtable->header, 1);
1913  ptyoff = muse_pfits_get_crval(aPixtable->header, 2);
1914  }
1915  cpl_propertylist_append_float(aPixtable->header, MUSE_HDR_PT_XLO, xlo + ptxoff);
1916  cpl_propertylist_append_float(aPixtable->header, MUSE_HDR_PT_XHI, xhi + ptxoff);
1917  cpl_propertylist_append_float(aPixtable->header, MUSE_HDR_PT_YLO, ylo + ptyoff);
1918  cpl_propertylist_append_float(aPixtable->header, MUSE_HDR_PT_YHI, yhi + ptyoff);
1919  cpl_propertylist_append_float(aPixtable->header, MUSE_HDR_PT_LLO, llo);
1920  cpl_propertylist_append_float(aPixtable->header, MUSE_HDR_PT_LHI, lhi);
1921  cpl_propertylist_append_int(aPixtable->header, MUSE_HDR_PT_ILO, ifulo);
1922  cpl_propertylist_append_int(aPixtable->header, MUSE_HDR_PT_IHI, ifuhi);
1923  cpl_propertylist_append_int(aPixtable->header, MUSE_HDR_PT_SLO, slicelo);
1924  cpl_propertylist_append_int(aPixtable->header, MUSE_HDR_PT_SHI, slicehi);
1925 
1926  return CPL_ERROR_NONE;
1927 } /* muse_pixtable_compute_limits() */
1928 
1929 /*---------------------------------------------------------------------------*/
1944 /*---------------------------------------------------------------------------*/
1945 cpl_error_code
1946 muse_pixtable_flux_multiply(muse_pixtable *aPixtable, double aScale)
1947 {
1948  cpl_ensure_code(aPixtable && aPixtable->table, CPL_ERROR_NULL_INPUT);
1949  cpl_errorstate prestate = cpl_errorstate_get();
1950  cpl_table_multiply_scalar(aPixtable->table, MUSE_PIXTABLE_DATA, aScale);
1951  cpl_table_multiply_scalar(aPixtable->table, MUSE_PIXTABLE_STAT, aScale*aScale);
1952  cpl_error_code rc = CPL_ERROR_NONE;
1953  if (!cpl_errorstate_is_equal(prestate)) {
1954  rc = cpl_error_get_code();
1955  }
1956  return rc;
1957 } /* muse_pixtable_flux_multiply() */
1958 
1959 /*---------------------------------------------------------------------------*/
1992 /*---------------------------------------------------------------------------*/
1993 cpl_error_code
1995  const cpl_array *aLambdas, const cpl_array *aFlux,
1996  muse_pixtable_operation aOperation)
1997 {
1998  cpl_ensure_code(aPixtable && aPixtable->table, CPL_ERROR_NULL_INPUT);
1999  cpl_ensure_code(aLambdas && aFlux, CPL_ERROR_NULL_INPUT);
2000  cpl_ensure_code(cpl_array_get_size(aLambdas) > 0 &&
2001  cpl_array_get_size(aLambdas) == cpl_array_get_size(aFlux),
2002  CPL_ERROR_ILLEGAL_INPUT);
2003  cpl_ensure_code(cpl_array_get_type(aLambdas) == CPL_TYPE_DOUBLE &&
2004  cpl_array_get_type(aFlux) == CPL_TYPE_DOUBLE,
2005  CPL_ERROR_INCOMPATIBLE_INPUT);
2006 
2007  /* To apply a given spectrum, decompose the full pixel table into per-slice *
2008  * pixel tables, which can be sorted by wavelength. For these, we can then *
2009  * linearly interpolate the input spectrum onto the same wavelengths as in *
2010  * the sub-pixel table and then subtract that in the form on an array. *
2011  * This also has the advantage of being parallelizable. */
2012  muse_pixtable **slicepts = muse_pixtable_extracted_get_slices(aPixtable);
2013  cpl_size islice, nslices = muse_pixtable_extracted_get_size(slicepts);
2014 
2015  switch (aOperation) {
2016  case MUSE_PIXTABLE_OPERATION_SUBTRACT:
2017  cpl_msg_debug(__func__, "Subtracting spectrum from pixel table with %"
2018  CPL_SIZE_FORMAT" slices...", nslices);
2019  break;
2020  case MUSE_PIXTABLE_OPERATION_MULTIPLY:
2021  cpl_msg_debug(__func__, "Multiplying pixel table of %"CPL_SIZE_FORMAT
2022  " slices with spectrum...", nslices);
2023  break;
2024  case MUSE_PIXTABLE_OPERATION_DIVIDE:
2025  cpl_msg_debug(__func__, "Dividing pixel table of %"CPL_SIZE_FORMAT
2026  " slices with spectrum...", nslices);
2027  break;
2028  default:
2030  return cpl_error_set(__func__, CPL_ERROR_UNSUPPORTED_MODE);
2031  } /* switch */
2032 
2033  /* now loop through all extracted slices in parallel */
2034  #pragma omp parallel for default(none) \
2035  shared(aFlux, aLambdas, aOperation, nslices, slicepts)
2036  for (islice = 0; islice < nslices; islice++) {
2037  muse_pixtable *thispt = slicepts[islice]; /* one slice pixel table */
2038  cpl_propertylist *order = cpl_propertylist_new();
2039  cpl_propertylist_append_bool(order, MUSE_PIXTABLE_LAMBDA, CPL_FALSE);
2040  cpl_table_sort(thispt->table, order);
2041  cpl_propertylist_delete(order);
2042 
2043  /* get sorted wavelengths as double array (so *
2044  * that the interpolation function can be used!) */
2045  cpl_table_cast_column(thispt->table, MUSE_PIXTABLE_LAMBDA, "lbda_d",
2046  CPL_TYPE_DOUBLE);
2047  cpl_array *lambdapt = muse_cpltable_extract_column(thispt->table, "lbda_d");
2048  cpl_array *fint = muse_cplarray_interpolate_linear(lambdapt, aLambdas,
2049  aFlux);
2050 
2051  /* get relevant pixel table columns to apply the flat-field spectrum to */
2052  cpl_array *datapt = muse_cpltable_extract_column(thispt->table,
2054  *statpt = muse_cpltable_extract_column(thispt->table,
2056  /* apply the flat-field spectrum */
2057  if (aOperation == MUSE_PIXTABLE_OPERATION_SUBTRACT) {
2058  cpl_array_subtract(datapt, fint);
2059  } else if (aOperation == MUSE_PIXTABLE_OPERATION_DIVIDE) {
2060  cpl_array_divide(datapt, fint);
2061  /* divide the variance squared (i.e. twice) */
2062  cpl_array_divide(statpt, fint);
2063  cpl_array_divide(statpt, fint);
2064  } else { /* MUSE_PIXTABLE_OPERATION_MULTIPLY */
2065  cpl_array_multiply(datapt, fint);
2066  /* multiply to the variance squared (i.e. twice) */
2067  cpl_array_multiply(statpt, fint);
2068  cpl_array_multiply(statpt, fint);
2069  } /* if aOperation */
2070  cpl_array_delete(fint);
2071 
2072  cpl_array_unwrap(datapt);
2073  cpl_array_unwrap(statpt);
2074  cpl_array_unwrap(lambdapt);
2075  cpl_table_erase_column(thispt->table, "lbda_d");
2076  } /* for islice (all slice pixel tables) */
2078 
2079  return CPL_ERROR_NONE;
2080 } /* muse_pixtable_spectrum_apply() */
2081 
2082 /*---------------------------------------------------------------------------*/
2092 /*---------------------------------------------------------------------------*/
2093 static cpl_error_code
2095 {
2096  cpl_ensure_code(aPixtable && aPixtable->header && aPixtable->table,
2097  CPL_ERROR_NULL_INPUT);
2098  if (cpl_table_count_selected(aPixtable->table) < 1) {
2099  return CPL_ERROR_NONE; /* nothing to do */
2100  }
2101  cpl_array *sel = cpl_table_where_selected(aPixtable->table);
2102  cpl_size narray = cpl_array_get_size(sel),
2103  nselprev = 0, /* selected ones in previous exposures */
2104  ifst = 0, ilst = 0;
2105  const cpl_size *asel = cpl_array_get_data_cplsize_const(sel);
2106  unsigned int nexp = 0;
2107  do {
2108  /* get exposure range for (next) exposure */
2109  char *kwfst = cpl_sprintf(MUSE_HDR_PT_EXP_FST, ++nexp),
2110  *kwlst = cpl_sprintf(MUSE_HDR_PT_EXP_LST, nexp);
2111  if (!cpl_propertylist_has(aPixtable->header, kwfst) ||
2112  !cpl_propertylist_has(aPixtable->header, kwlst)) {
2113  cpl_free(kwfst);
2114  cpl_free(kwlst);
2115  break; /* nothing (more) to do */
2116  }
2117  ifst = cpl_propertylist_get_long_long(aPixtable->header, kwfst);
2118  ilst = cpl_propertylist_get_long_long(aPixtable->header, kwlst);
2119  /* count selected */
2120  cpl_size i, nsel = 0;
2121  for (i = 0; i < narray; i++) {
2122  if (asel[i] >= ifst && asel[i] <= ilst) {
2123  nsel++;
2124  } /* if in range of this exposure */
2125  } /* for i (all array entries) */
2126  cpl_size ifst2 = ifst - nselprev,
2127  ilst2 = ilst - nsel - nselprev;
2128  cpl_msg_debug(__func__, "exp %d old %"CPL_SIZE_FORMAT"..%"CPL_SIZE_FORMAT
2129  ", %"CPL_SIZE_FORMAT" selected (previous: %"CPL_SIZE_FORMAT
2130  "), new %"CPL_SIZE_FORMAT"..%"CPL_SIZE_FORMAT, nexp,
2131  ifst, ilst, nsel, nselprev, ifst2, ilst2);
2132  /* finally update the header entries */
2133  muse_cplpropertylist_update_long_long(aPixtable->header, kwfst, ifst2);
2134  muse_cplpropertylist_update_long_long(aPixtable->header, kwlst, ilst2);
2135  cpl_free(kwfst);
2136  cpl_free(kwlst);
2137  nselprev += nsel; /* add selected ones for previous exposures */
2138  } while (ilst >= ifst);
2139  cpl_array_delete(sel);
2140  return CPL_ERROR_NONE;
2141 } /* muse_pixtable_fix_exp_headers() */
2142 
2143 /*---------------------------------------------------------------------------*/
2161 /*---------------------------------------------------------------------------*/
2162 cpl_error_code
2164  double aHigh)
2165 {
2166  cpl_ensure_code(aPixtable && aPixtable->table && aPixtable->header,
2167  CPL_ERROR_NULL_INPUT);
2168  if (cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LLO) > aLow &&
2169  cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LHI) < aHigh) {
2170  /* no need to do anything */
2171  return CPL_ERROR_NONE;
2172  }
2173 #pragma omp critical(cpl_table_select)
2174  {
2175  cpl_table_unselect_all(aPixtable->table);
2176  cpl_table_or_selected_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA,
2177  CPL_LESS_THAN, aLow);
2178  cpl_table_or_selected_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA,
2179  CPL_GREATER_THAN, aHigh);
2180  muse_pixtable_fix_exp_headers(aPixtable);
2181  cpl_table_erase_selected(aPixtable->table);
2182  }
2183 
2184  /* do the same for the flat-field spectrum, if it exists, but leave a *
2185  * spectral resolution element more on each side, for interpolation */
2186 #pragma omp critical(cpl_table_select)
2187  if (aPixtable->ffspec) {
2188  cpl_table_unselect_all(aPixtable->ffspec);
2189  cpl_table_or_selected_double(aPixtable->ffspec, MUSE_PIXTABLE_FFLAMBDA,
2190  CPL_LESS_THAN,
2191  aLow - 2. * kMuseSpectralSamplingA);
2192  cpl_table_or_selected_double(aPixtable->ffspec, MUSE_PIXTABLE_FFLAMBDA,
2193  CPL_GREATER_THAN,
2194  aHigh + 2. * kMuseSpectralSamplingA);
2195  cpl_table_erase_selected(aPixtable->ffspec);
2196  }
2197 
2198  return muse_pixtable_compute_limits(aPixtable);
2199 } /* muse_pixtable_restrict_wavelength() */
2200 
2201 /*---------------------------------------------------------------------------*/
2212 /*---------------------------------------------------------------------------*/
2213 cpl_error_code
2214 muse_pixtable_restrict_xpos(muse_pixtable *aPixtable, double aLo, double aHi)
2215 {
2216  cpl_ensure_code(aPixtable && aPixtable->table && aPixtable->header,
2217  CPL_ERROR_NULL_INPUT);
2218  if (cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_XLO) > aLo &&
2219  cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_XHI) < aHi) {
2220  /* no need to do anything */
2221  return CPL_ERROR_NONE;
2222  }
2223  double ptxoff = 0.;
2225  /* need to use the real coordinate offset for celestial spherical */
2226  ptxoff = muse_pfits_get_crval(aPixtable->header, 1);
2227  }
2228 #pragma omp critical(cpl_table_select)
2229  {
2230  cpl_table_unselect_all(aPixtable->table);
2231  cpl_table_or_selected_float(aPixtable->table, MUSE_PIXTABLE_XPOS,
2232  CPL_LESS_THAN, aLo - ptxoff);
2233  cpl_table_or_selected_float(aPixtable->table, MUSE_PIXTABLE_XPOS,
2234  CPL_GREATER_THAN, aHi - ptxoff);
2235  muse_pixtable_fix_exp_headers(aPixtable);
2236  cpl_table_erase_selected(aPixtable->table);
2237  }
2238  return muse_pixtable_compute_limits(aPixtable);
2239 } /* muse_pixtable_restrict_xpos() */
2240 
2241 /*---------------------------------------------------------------------------*/
2252 /*---------------------------------------------------------------------------*/
2253 cpl_error_code
2254 muse_pixtable_restrict_ypos(muse_pixtable *aPixtable, double aLo, double aHi)
2255 {
2256  cpl_ensure_code(aPixtable && aPixtable->table && aPixtable->header,
2257  CPL_ERROR_NULL_INPUT);
2258  if (cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_YLO) > aLo &&
2259  cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_YHI) < aHi) {
2260  /* no need to do anything */
2261  return CPL_ERROR_NONE;
2262  }
2263  double ptyoff = 0.;
2265  /* need to use the real coordinate offset for celestial spherical */
2266  ptyoff = muse_pfits_get_crval(aPixtable->header, 2);
2267  }
2268 #pragma omp critical(cpl_table_select)
2269  {
2270  cpl_table_unselect_all(aPixtable->table);
2271  cpl_table_or_selected_float(aPixtable->table, MUSE_PIXTABLE_YPOS,
2272  CPL_LESS_THAN, aLo - ptyoff);
2273  cpl_table_or_selected_float(aPixtable->table, MUSE_PIXTABLE_YPOS,
2274  CPL_GREATER_THAN, aHi - ptyoff);
2275  muse_pixtable_fix_exp_headers(aPixtable);
2276  cpl_table_erase_selected(aPixtable->table);
2277  }
2278  return muse_pixtable_compute_limits(aPixtable);
2279 } /* muse_pixtable_restrict_ypos() */
2280 
2281 /*---------------------------------------------------------------------------*/
2295 /*---------------------------------------------------------------------------*/
2296 cpl_error_code
2297 muse_pixtable_erase_ifu_slice(muse_pixtable *aPixtable, unsigned char aIFU,
2298  unsigned short aSlice)
2299 {
2300  cpl_ensure_code(aPixtable, CPL_ERROR_NULL_INPUT);
2301  cpl_size nrow = muse_pixtable_get_nrow(aPixtable);
2302  cpl_ensure_code(nrow > 0, CPL_ERROR_DATA_NOT_FOUND);
2303 
2304  cpl_table_unselect_all(aPixtable->table);
2305  uint32_t *origin = (uint32_t *)cpl_table_get_data_int(aPixtable->table,
2307  cpl_size irow;
2308  for (irow = 0; irow < nrow; irow++) {
2309  unsigned char ifu = muse_pixtable_origin_get_ifu(origin[irow]);
2310  unsigned short slice = muse_pixtable_origin_get_slice(origin[irow]);
2311  if (ifu == aIFU && slice == aSlice) {
2312  cpl_table_select_row(aPixtable->table, irow);
2313  } /* if same IFU and slice */
2314  } /* for irow (all pixtable rows) */
2315  cpl_size nsel = cpl_table_count_selected(aPixtable->table);
2316  cpl_error_code rc = cpl_table_erase_selected(aPixtable->table);
2317  cpl_msg_debug(__func__, "Erased %"CPL_SIZE_FORMAT" rows from pixel table",
2318  nsel);
2319 
2320  muse_pixtable_fix_exp_headers(aPixtable);
2321  muse_pixtable_compute_limits(aPixtable);
2322  return rc;
2323 } /* muse_pixtable_erase_ifu_slice() */
2324 
2325 /*---------------------------------------------------------------------------*/
2342 /*---------------------------------------------------------------------------*/
2343 cpl_error_code
2345 {
2346  cpl_ensure_code(aPixtable && aPixtable->table, CPL_ERROR_NULL_INPUT);
2347  cpl_ensure_code(aMask && aMask->mask, CPL_ERROR_NULL_INPUT);
2348  float *cdata_xpos = cpl_table_get_data_float(aPixtable->table,
2350  float *cdata_ypos = cpl_table_get_data_float(aPixtable->table,
2352  cpl_size n_rows = cpl_table_get_nrow(aPixtable->table);
2353  cpl_size i_row;
2354  double crval_x = 0;
2355  double crpix_x = 1;
2356  double cdelt_x = 1;
2357  double crval_y = 0;
2358  double crpix_y = 1;
2359  double cdelt_y = 1;
2360  if (aMask->header) {
2361  crval_x = muse_pfits_get_crval(aMask->header, 1);
2362  crpix_x = muse_pfits_get_crpix(aMask->header, 1);
2363  cdelt_x = muse_pfits_get_cd(aMask->header, 1, 1);
2364  crval_y = muse_pfits_get_crval(aMask->header, 2);
2365  crpix_y = muse_pfits_get_crpix(aMask->header, 2);
2366  cdelt_y = muse_pfits_get_cd(aMask->header, 2, 2);
2367  }
2368  cpl_size nx = cpl_mask_get_size_x(aMask->mask);
2369  cpl_size ny = cpl_mask_get_size_y(aMask->mask);
2370  cpl_size n_enabled = cpl_mask_count(aMask->mask);
2371  cpl_msg_debug(__func__, "Mask contains %"CPL_SIZE_FORMAT" (%.2f %%) enabled "
2372  "pixels of %"CPL_SIZE_FORMAT" total", n_enabled,
2373  100.*n_enabled/nx/ny, nx*ny);
2374  cpl_size n_sel = n_rows;
2375  cpl_size n_in_table = 0;
2376  for (i_row = 0; i_row < n_rows; i_row++) {
2377  cpl_size ix = lround((cdata_xpos[i_row] - crval_x) / cdelt_x + crpix_x);
2378  cpl_size iy = lround((cdata_ypos[i_row] - crval_y) / cdelt_y + crpix_y);
2379  if ((ix < 1) || (ix > nx) || (iy < 1) || (iy > ny)) {
2380  continue;
2381  }
2382  n_in_table++;
2383  if (cpl_mask_get(aMask->mask, ix, iy) != CPL_BINARY_1) {
2384  cpl_table_unselect_row(aPixtable->table, i_row);
2385  n_sel--;
2386  }
2387  }
2388  cpl_msg_debug(__func__, "Mask selected %"CPL_SIZE_FORMAT" (%.2f %%/%.2f %%) "
2389  "pixels of %"CPL_SIZE_FORMAT" total/%"CPL_SIZE_FORMAT" in mask "
2390  "area", n_sel, 100.*n_sel/n_rows, 100.*n_sel/n_in_table,
2391  n_rows, n_in_table);
2392  return CPL_ERROR_NONE;
2393 } /* muse_pixtable_and_selected_mask() */
2394 
2395 /*---------------------------------------------------------------------------*/
2422 /*---------------------------------------------------------------------------*/
2423 cpl_error_code
2424 muse_pixtable_dump(muse_pixtable *aPixtable, cpl_size aStart, cpl_size aCount,
2425  unsigned char aDisplayHeader)
2426 {
2427  cpl_ensure_code(aPixtable && aPixtable->table && aPixtable->header,
2428  CPL_ERROR_NULL_INPUT);
2429  cpl_size nrows = muse_pixtable_get_nrow(aPixtable);
2430  cpl_ensure_code(aStart >= 0 && aStart < nrows && aCount >= 0,
2431  CPL_ERROR_ILLEGAL_INPUT);
2432  cpl_size last = aStart + aCount;
2433  if (last > nrows - 1) {
2434  last = nrows;
2435  }
2436  int haswcs = muse_pixtable_wcs_check(aPixtable);
2437  double ptxoff = 0., ptyoff = 0.;
2438  if (haswcs == MUSE_PIXTABLE_WCS_CELSPH) {
2439  /* need to use the real coordinate offset for celestial spherical */
2440  ptxoff = muse_pfits_get_crval(aPixtable->header, 1);
2441  ptyoff = muse_pfits_get_crval(aPixtable->header, 2);
2442  }
2443  /* for non-pixel coordinates, we will need different formats below */
2444  haswcs = (haswcs == MUSE_PIXTABLE_WCS_NATSPH ||
2445  haswcs == MUSE_PIXTABLE_WCS_CELSPH);
2446  float *cdata_xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
2447  *cdata_ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS),
2448  *cdata_lambda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA),
2449  *cdata_data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA),
2450  *cdata_stat = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_STAT);
2451  cpl_errorstate es = cpl_errorstate_get();
2452  float *cdata_weight = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_WEIGHT);
2453  cpl_errorstate_set(es); /* ignore errors due to missing weight column */
2454  int *cdata_dq = cpl_table_get_data_int(aPixtable->table, MUSE_PIXTABLE_DQ);
2455  uint32_t *cdata_origin = (uint32_t *)cpl_table_get_data_int(aPixtable->table,
2457  cpl_ensure_code(cdata_xpos && cdata_ypos && cdata_lambda &&
2458  cdata_data && cdata_dq && cdata_stat,
2459  CPL_ERROR_BAD_FILE_FORMAT);
2460 
2461  /* print header */
2462  if (aDisplayHeader) {
2463  printf("# xpos ypos lambda data dq stat"
2464  " weight exposure IFU xCCD yCCD xRaw yRaw slice\n");
2465  }
2466  if (aDisplayHeader == 1) {
2467  printf("#%13s %13s %9s %11s flag %11s ---------- No No pix "
2468  " pix pix pix No\n# flux in [%s]\n# flux**2 in [%s]\n",
2469  cpl_table_get_column_unit(aPixtable->table, MUSE_PIXTABLE_XPOS),
2470  cpl_table_get_column_unit(aPixtable->table, MUSE_PIXTABLE_YPOS),
2471  cpl_table_get_column_unit(aPixtable->table, MUSE_PIXTABLE_LAMBDA),
2472  "(flux)", "(flux**2)",
2473  cpl_table_get_column_unit(aPixtable->table, MUSE_PIXTABLE_DATA),
2474  cpl_table_get_column_unit(aPixtable->table, MUSE_PIXTABLE_STAT));
2475  }
2476 
2477  cpl_size i;
2478  for (i = aStart; i < last; i++) {
2479  int x = muse_pixtable_origin_get_x(cdata_origin[i], aPixtable, i),
2480  y = muse_pixtable_origin_get_y_fast(cdata_origin[i]),
2481  xraw = x, yraw = y;
2482  muse_quadrants_coords_to_raw(NULL, &xraw, &yraw);
2483  if (haswcs) {
2484  printf("%14.7e %14.7e %9.3f ", cdata_xpos[i] + ptxoff, cdata_ypos[i] + ptyoff,
2485  cdata_lambda[i]);
2486  } else {
2487  printf("%14.8f %14.8f %9.3f ", cdata_xpos[i], cdata_ypos[i], cdata_lambda[i]);
2488  }
2489  printf("%12.5e 0x%08x %11.5e %10.4e %2d %2d %4d %4d %4d %4d %2d\n",
2490  cdata_data[i], cdata_dq[i], cdata_stat[i],
2491  cdata_weight ? cdata_weight[i] : 0.,
2492  muse_pixtable_get_expnum(aPixtable, i),
2493  cdata_origin ? muse_pixtable_origin_get_ifu_fast(cdata_origin[i])
2494  : 0,
2495  x, y, xraw, yraw,
2496  cdata_origin ? muse_pixtable_origin_get_slice_fast(cdata_origin[i])
2497  : 0);
2498  } /* for i (pixel table rows) */
2499 
2500  return CPL_ERROR_NONE;
2501 } /* muse_pixtable_dump() */
2502 
2503 /*---------------------------------------------------------------------------*/
2525 /*---------------------------------------------------------------------------*/
2528 {
2529  cpl_ensure(aPixtable, CPL_ERROR_NULL_INPUT, MUSE_PIXTABLE_WCS_UNKNOWN);
2530  const char *unitx = cpl_table_get_column_unit(aPixtable->table,
2532  *unity = cpl_table_get_column_unit(aPixtable->table,
2534  cpl_ensure(unitx, CPL_ERROR_DATA_NOT_FOUND, MUSE_PIXTABLE_WCS_UNKNOWN);
2535  cpl_ensure(unity, CPL_ERROR_DATA_NOT_FOUND, MUSE_PIXTABLE_WCS_UNKNOWN);
2536  /* should be equal in the first three characters */
2537  cpl_ensure(!strncmp(unitx, unity, 4), CPL_ERROR_INCOMPATIBLE_INPUT,
2539  if (!strncmp(unitx, "deg", 4)) {
2540  return MUSE_PIXTABLE_WCS_CELSPH;
2541  }
2542  if (!strncmp(unitx, "pix", 4)) {
2543  return MUSE_PIXTABLE_WCS_PIXEL;
2544  }
2545  if (!strncmp(unitx, "rad", 4)) {
2546  return MUSE_PIXTABLE_WCS_NATSPH;
2547  }
2548  cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
2550 } /* muse_pixtable_wcs_check() */
2551 
2552 /*---------------------------------------------------------------------------*/
2560 /*---------------------------------------------------------------------------*/
2561 cpl_boolean
2563 {
2564  cpl_ensure(aPixtable, CPL_ERROR_NULL_INPUT, CPL_FALSE);
2565  cpl_errorstate prestate = cpl_errorstate_get();
2566  cpl_boolean flag = cpl_propertylist_get_bool(aPixtable->header,
2568  cpl_errorstate_set(prestate);
2569  return flag;
2570 }
2571 
2572 /*---------------------------------------------------------------------------*/
2580 /*---------------------------------------------------------------------------*/
2581 cpl_boolean
2583 {
2584  cpl_ensure(aPixtable, CPL_ERROR_NULL_INPUT, CPL_FALSE);
2585  cpl_errorstate prestate = cpl_errorstate_get();
2586  cpl_boolean flag = cpl_propertylist_get_bool(aPixtable->header,
2588  cpl_errorstate_set(prestate);
2589  return flag;
2590 }
2591 
2592 /*---------------------------------------------------------------------------*/
2600 /*---------------------------------------------------------------------------*/
2601 cpl_boolean
2603 {
2604  cpl_ensure(aPixtable, CPL_ERROR_NULL_INPUT, CPL_FALSE);
2605  return cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_RVCORR);
2606 } /* muse_pixtable_is_rvcorr() */
2607 
2608 /*---------------------------------------------------------------------------*/
2617 /*---------------------------------------------------------------------------*/
2618 cpl_error_code
2619 muse_pixtable_reset_dq(muse_pixtable *aPixtable, unsigned int aDQ)
2620 {
2621  cpl_ensure_code(aPixtable, CPL_ERROR_NULL_INPUT);
2622 
2623  unsigned int *dq = (unsigned int *)cpl_table_get_data_int(aPixtable->table,
2625  inverse = ~aDQ; /* inverse bitmask to AND each row's DQ with */
2626  cpl_size i, nrow = muse_pixtable_get_nrow(aPixtable);
2627  #pragma omp parallel for default(none) /* as req. by Ralf */ \
2628  shared(dq, inverse, nrow)
2629  for (i = 0; i < nrow; i++) {
2630  dq[i] &= inverse;
2631  } /* for i (all table rows) */
2632  return CPL_ERROR_NONE;
2633 } /* muse_pixtable_reset_dq() */
2634 
2635 /*---------------------------------------------------------------------------*/
2652 /*---------------------------------------------------------------------------*/
2655 {
2656  cpl_ensure(aPixtable && aPixtable->header, CPL_ERROR_NULL_INPUT, NULL);
2657  unsigned int expnum = muse_pixtable_get_expnum(aPixtable, 0),
2658  explast = muse_pixtable_get_expnum(aPixtable,
2659  muse_pixtable_get_nrow(aPixtable) - 1);
2660  cpl_ensure(expnum == explast, CPL_ERROR_ILLEGAL_INPUT, NULL);
2661 
2662  /* data seems to be valid, we can create the output list */
2664 
2665  /* split the pixel table up into per-slice pixel tables to handle *
2666  * them separately (to get the x-range of each slice for the header) */
2668  /* variables for current image and current IFU */
2669  muse_image *image = NULL;
2670  unsigned short ifu = 0, /* we haven't found any IFU yet */
2671  ilist = 0; /* image index in the list */
2672  int ipt, npt = muse_pixtable_extracted_get_size(pts);
2673  for (ipt = 0; ipt < npt; ipt++) {
2674  float *cdata = cpl_table_get_data_float(pts[ipt]->table, MUSE_PIXTABLE_DATA),
2675  *cstat = cpl_table_get_data_float(pts[ipt]->table, MUSE_PIXTABLE_STAT);
2676  int *cdq = cpl_table_get_data_int(pts[ipt]->table, MUSE_PIXTABLE_DQ);
2677  uint32_t *corigin = (uint32_t *)cpl_table_get_data_int(pts[ipt]->table,
2679  /* if we got to the next (or first) IFU, create new output *
2680  * image of the size of a typical MUSE CCD copy the header *
2681  * from the pixel table but remove specific entries */
2682  if (ifu != muse_pixtable_origin_get_ifu_fast(corigin[0])) {
2683  image = muse_image_new();
2684  image->header = cpl_propertylist_duplicate(pts[ipt]->header);
2685  cpl_propertylist_erase_regexp(image->header, "^ESO DRS MUSE PIXTABLE", 0);
2686  image->data = cpl_image_new(kMuseOutputXRight, kMuseOutputYTop, CPL_TYPE_FLOAT);
2687  image->dq = cpl_image_new(kMuseOutputXRight, kMuseOutputYTop, CPL_TYPE_INT);
2688  /* fill DQ image with EURO3D_MISSDATA, upper value will be cast correctly */
2689  cpl_image_fill_noise_uniform(image->dq, EURO3D_MISSDATA, EURO3D_MISSDATA + 0.1);
2690  image->stat = cpl_image_new(kMuseOutputXRight, kMuseOutputYTop, CPL_TYPE_FLOAT);
2691  cpl_msg_debug(__func__, "new image (index %hu in list)", ilist);
2692  muse_imagelist_set(list, image, ilist++);
2693  } /* if ifu */
2694  if (!image) { /* it cannot really go wrong, but to be sure... */
2695  cpl_msg_error(__func__, "ipt = %d: no image!", ipt);
2696  continue;
2697  }
2698  float *idata = cpl_image_get_data_float(image->data),
2699  *istat = cpl_image_get_data_float(image->stat);
2700  int *idq = cpl_image_get_data_int(image->dq);
2701 
2702  ifu = muse_pixtable_origin_get_ifu_fast(corigin[0]);
2703  unsigned short slice = muse_pixtable_origin_get_slice_fast(corigin[0]);
2704  unsigned int xoff = muse_pixtable_origin_get_offset(pts[ipt], expnum, ifu,
2705  slice),
2706  x1 = INT_MAX, x2 = 0,
2707  irow, nrow = muse_pixtable_get_nrow(pts[ipt]);
2708  for (irow = 0; irow < nrow; irow++) {
2709  /* get coordinate indices from the origin column */
2710  unsigned int x = muse_pixtable_origin_get_x_fast(corigin[irow], xoff) - 1,
2711  y = muse_pixtable_origin_get_y_fast(corigin[irow]) - 1;
2712  idata[x + y*kMuseOutputXRight] = cdata[irow];
2713  idq[x + y*kMuseOutputXRight] = cdq[irow];
2714  istat[x + y*kMuseOutputXRight] = cstat[irow];
2715  if (x < x1) {
2716  x1 = x;
2717  }
2718  if (x > x2) {
2719  x2 = x;
2720  }
2721  } /* for irow */
2722  /* record the approximate horizontal center in the output header */
2723  char *keyword = cpl_sprintf("ESO DRS MUSE SLICE%hu CENTER", slice);
2724  cpl_propertylist_update_float(image->header, keyword, (x2 + x1) / 2. + 1.);
2725 #if 0
2726  cpl_msg_debug(__func__, "IFU %hu %s = %.1f", ifu, keyword,
2727  (x2 + x1) / 2. + 1.);
2728 #endif
2729  cpl_free(keyword);
2730  } /* for ipt */
2732 
2733  return list;
2734 } /* muse_pixtable_to_imagelist() */
2735 
2736 /*---------------------------------------------------------------------------*/
2762 /*---------------------------------------------------------------------------*/
2763 cpl_error_code
2765 {
2766  cpl_ensure_code(aPixtable && aPixtable->header && aList, CPL_ERROR_NULL_INPUT);
2767  unsigned int expnum = muse_pixtable_get_expnum(aPixtable, 0),
2768  explast = muse_pixtable_get_expnum(aPixtable,
2769  muse_pixtable_get_nrow(aPixtable) - 1);
2770  cpl_ensure_code(expnum == explast, CPL_ERROR_ILLEGAL_INPUT);
2771 
2772  /* split the pixel table up into per-slice pixel tables *
2773  * as in muse_pixtable_to_imagelist() and ensure that the *
2774  * number of IFUs matches the number of input images */
2776  if (muse_pixtable_extracted_get_size(pts) / kMuseSlicesPerCCD
2777  != muse_imagelist_get_size(aList)) {
2779  return cpl_error_set(__func__, CPL_ERROR_INCOMPATIBLE_INPUT);
2780  }
2781  /* variables for current image and current IFU */
2782  muse_image *image = NULL;
2783  unsigned short ifu = 0, /* we haven't worked with any IFU yet */
2784  ilist = 0; /* image index in the list */
2785  int ipt, npt = muse_pixtable_extracted_get_size(pts);
2786  for (ipt = 0; ipt < npt; ipt++) {
2787  float *cdata = cpl_table_get_data_float(pts[ipt]->table, MUSE_PIXTABLE_DATA),
2788  *cstat = cpl_table_get_data_float(pts[ipt]->table, MUSE_PIXTABLE_STAT);
2789  uint32_t *corigin = (uint32_t *)cpl_table_get_data_int(pts[ipt]->table,
2791  /* if we got to the next (or first) IFU, get an image from the list */
2792  if (ifu != muse_pixtable_origin_get_ifu_fast(corigin[0])) {
2793  image = muse_imagelist_get(aList, ilist++);
2794  } /* if ifu */
2795  if (!image) { /* it cannot really go wrong, but to be sure... */
2796  cpl_msg_error(__func__, "ipt = %d: no image!", ipt);
2797  continue;
2798  }
2799  float *idata = cpl_image_get_data_float(image->data),
2800  *istat = cpl_image_get_data_float(image->stat);
2801  ifu = muse_pixtable_origin_get_ifu_fast(corigin[0]);
2802  unsigned short slice = muse_pixtable_origin_get_slice_fast(corigin[0]);
2803  unsigned int xoff = muse_pixtable_origin_get_offset(pts[ipt], expnum, ifu,
2804  slice),
2805  irow, nrow = muse_pixtable_get_nrow(pts[ipt]);
2806  for (irow = 0; irow < nrow; irow++) {
2807  /* get coordinate indices from the origin column */
2808  unsigned int x = muse_pixtable_origin_get_x_fast(corigin[irow], xoff) - 1,
2809  y = muse_pixtable_origin_get_y_fast(corigin[irow]) - 1;
2810  cdata[irow] = idata[x + y*kMuseOutputXRight];
2811  cstat[irow] = istat[x + y*kMuseOutputXRight];
2812  } /* for irow */
2813  } /* for ipt */
2815 
2816  return CPL_ERROR_NONE;
2817 } /* muse_pixtable_from_imagelist() */
2818 
2819 /*---------------------------------------------------------------------------*/
2832 /*---------------------------------------------------------------------------*/
2833 muse_pixtable **
2835 {
2836  cpl_ensure(aPixtable, CPL_ERROR_NULL_INPUT, NULL);
2837  cpl_size n_rows = cpl_table_get_nrow(aPixtable->table);
2838  unsigned int ifu_slice_mask = (0x1f << MUSE_ORIGIN_SHIFT_IFU) | 0x3f;
2839  cpl_table_duplicate_column(aPixtable->table, "ifuslice",
2840  aPixtable->table, MUSE_PIXTABLE_ORIGIN);
2841  unsigned int *slicedata = (unsigned int *)
2842  cpl_table_get_data_int(aPixtable->table, "ifuslice");
2843  cpl_size i_row;
2844  unsigned int last_ifu_slice = 0;
2845  int is_sorted = CPL_TRUE;
2846  for (i_row = 0; i_row < n_rows; i_row++) {
2847  slicedata[i_row] &= ifu_slice_mask;
2848  if (is_sorted && slicedata[i_row] < last_ifu_slice) {
2849  is_sorted = CPL_FALSE;
2850  } else {
2851  last_ifu_slice = slicedata[i_row];
2852  }
2853  }
2854  if (!is_sorted) {
2855  cpl_propertylist *order = cpl_propertylist_new();
2856  cpl_propertylist_append_bool(order, "ifuslice", CPL_FALSE);
2857  cpl_propertylist_append_bool(order, MUSE_PIXTABLE_LAMBDA, CPL_FALSE);
2858  cpl_msg_debug(__func__, "sorting pixel table: quick sort, %"CPL_SIZE_FORMAT
2859  " entries", n_rows);
2860  cpl_table_sort(aPixtable->table, order);
2861  cpl_propertylist_delete(order);
2862  /* erase headers that depend on the order in the pixel table */
2863  cpl_propertylist_erase_regexp(aPixtable->header, MUSE_HDR_PT_EXP_REGEXP, 0);
2864  cpl_msg_debug(__func__, "pixel table sorted.");
2865  } /* if !sorted */
2866 
2867  i_row = 0;
2868  cpl_size n_col = cpl_table_get_ncol(aPixtable->table);
2869  cpl_array *colnames = cpl_table_get_column_names(aPixtable->table);
2870  muse_pixtable **slice_tables = cpl_calloc(1, sizeof(muse_pixtable *));
2871  cpl_size n_slices = 0;
2872  while (i_row < n_rows) {
2873  unsigned int ifu_slice = slicedata[i_row];
2874  cpl_size j_row;
2875  for (j_row = i_row+1; j_row < n_rows && slicedata[j_row] == ifu_slice;
2876  j_row++)
2877  ;
2878  cpl_size nrows_slice = j_row - i_row;
2879  muse_pixtable *slice_pixtable = cpl_calloc(1, sizeof(muse_pixtable));
2880  slice_pixtable->table = cpl_table_new(nrows_slice);
2881  cpl_size i_col;
2882  for (i_col = 0; i_col < n_col; i_col++) {
2883  const char *cname = cpl_array_get_string(colnames, i_col);
2884  if (strcmp(cname, "ifuslice") == 0)
2885  continue;
2886  cpl_type ctype = cpl_table_get_column_type(aPixtable->table, cname);
2887  if (ctype == CPL_TYPE_INT) {
2888  int *cdata = cpl_table_get_data_int(aPixtable->table, cname);
2889  cpl_table_wrap_int(slice_pixtable->table, cdata + i_row, cname);
2890  } else if (ctype == CPL_TYPE_FLOAT) {
2891  float *cdata = cpl_table_get_data_float(aPixtable->table, cname);
2892  cpl_table_wrap_float(slice_pixtable->table, cdata + i_row, cname);
2893  } else if (ctype == CPL_TYPE_DOUBLE) {
2894  double *cdata = cpl_table_get_data_double(aPixtable->table, cname);
2895  cpl_table_wrap_double(slice_pixtable->table, cdata + i_row, cname);
2896  } else if (ctype == CPL_TYPE_STRING) {
2897  char **cdata = cpl_table_get_data_string(aPixtable->table, cname);
2898  cpl_table_wrap_string(slice_pixtable->table, cdata + i_row, cname);
2899  }
2900  const char *unit = cpl_table_get_column_unit(aPixtable->table, cname);
2901  cpl_table_set_column_unit(slice_pixtable->table, cname, unit);
2902  } /* for i_col (all table columns) */
2903 
2904  slice_pixtable->header = cpl_propertylist_duplicate(aPixtable->header);
2905  muse_pixtable_compute_limits(slice_pixtable);
2906  slice_tables = cpl_realloc(slice_tables,
2907  (n_slices + 2) * sizeof(muse_pixtable *));
2908  slice_tables[n_slices] = slice_pixtable;
2909  n_slices++;
2910  slice_tables[n_slices] = NULL;
2911  i_row = j_row;
2912  } /* while (all full pixel table rows) */
2913  cpl_array_delete(colnames);
2914  cpl_table_erase_column(aPixtable->table, "ifuslice");
2915 
2916  return slice_tables;
2917 } /* muse_pixtable_extracted_get_slices() */
2918 
2919 /*---------------------------------------------------------------------------*/
2930 /*---------------------------------------------------------------------------*/
2931 cpl_size
2933 {
2934  cpl_ensure(aPixtables, CPL_ERROR_NULL_INPUT, -1);
2935  cpl_size n = 0;
2936  while (aPixtables[n] != NULL) {
2937  n++;
2938  }
2939  return n;
2940 } /* muse_pixtable_extracted_get_size() */
2941 
2942 /*---------------------------------------------------------------------------*/
2952 /*---------------------------------------------------------------------------*/
2953 void
2955 {
2956  if (!aPixtables) {
2957  return;
2958  }
2959  muse_pixtable **t;
2960  for (t = aPixtables; *t != NULL; t++) {
2961  cpl_array *colnames = cpl_table_get_column_names((*t)->table);
2962  cpl_size n_col = cpl_table_get_ncol((*t)->table);
2963  cpl_size i_col;
2964  for (i_col = 0; i_col < n_col; i_col++) {
2965  const char *cname = cpl_array_get_string(colnames, i_col);
2966  cpl_table_unwrap((*t)->table, cname);
2967  }
2968  cpl_array_delete(colnames);
2969  cpl_table_delete((*t)->table);
2970  cpl_propertylist_delete((*t)->header);
2971  cpl_free(*t);
2972  }
2973  cpl_free(aPixtables);
2974 } /* muse_pixtable_extracted_delete() */
2975 
#define MUSE_PIXTABLE_DQ
Definition: muse_pixtable.h:49
#define MUSE_PIXTABLE_XPOS
Definition: muse_pixtable.h:51
cpl_polynomial ** muse_trace_table_get_polys_for_slice(const cpl_table *aTable, const unsigned short aSlice)
construct polynomial from the trace table entry for the given slice
#define MUSE_PIXTABLE_FFDATA
Definition: muse_pixtable.h:61
muse_pixtable_wcs
State of the astrometric calibration of a MUSE pixel table.
unsigned int muse_pixtable_get_expnum(muse_pixtable *aPixtable, cpl_size aRow)
Get the exposure number of a given row in a pixel table.
cpl_array * muse_cplarray_interpolate_table_linear(const cpl_array *aTargetAbscissa, const cpl_table *aSrcTable, const char *aSrcAbscissa, const char *aSrcOrdinate)
Linear interpolation of a 1d column.
cpl_error_code muse_wave_table_get_orders(const cpl_table *aWave, unsigned short *aXOrder, unsigned short *aYOrder)
Determine the x- and y-order of the polynomial stored in a wavelength calibration table...
Structure definition for a collection of muse_images.
int muse_trace_table_get_order(const cpl_table *aTable)
determine order of tracing polynomial from table
#define MUSE_HDR_PT_EXP_FST
FITS header keyword defining the first row index for a given exposure.
Definition: muse_pixtable.h:93
void muse_pixtable_extracted_delete(muse_pixtable **aPixtables)
Delete a pixel table array.
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
cpl_polynomial * muse_wave_table_get_poly_for_slice(const cpl_table *aTable, const unsigned short aSlice)
Construct polynomial from the wavelength calibration table entry for the given slice.
const char * muse_pfits_get_extname(const cpl_propertylist *aHeaders)
find out the extension name
Definition: muse_pfits.c:205
muse_pixtable * muse_pixtable_load(const char *aFilename)
Load the table itself and the FITS headers of a MUSE pixel table from a file.
unsigned short muse_pixtable_origin_get_slice(uint32_t aOrigin)
Get the slice number from the encoded 32bit origin number.
muse_pixtable * muse_pixtable_duplicate(muse_pixtable *aPixtable)
Make a copy of the pixtanle.
#define MUSE_HDR_PT_XLO
FITS header keyword contains the lower limit of the data in x-direction.
cpl_size muse_pixtable_extracted_get_size(muse_pixtable **aPixtables)
Get the size of an array of extracted pixel tables.
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
unsigned int muse_pixtable_origin_set_offset(muse_pixtable *aPixtable, cpl_polynomial *aLTrace, unsigned short aIFU, unsigned short aSlice)
Set the slice offset from the pixel table header.
unsigned char muse_utils_get_ifu(const cpl_propertylist *aHeaders)
Find out the IFU/channel from which this header originated.
Definition: muse_utils.c:98
cpl_size muse_pfits_get_naxis(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the size of a given axis
Definition: muse_pfits.c:242
#define MUSE_HDR_PT_FLUXCAL
muse_pixtable ** muse_pixtable_extracted_get_slices(muse_pixtable *aPixtable)
Extract one pixel table per IFU and slice.
cpl_image * data
the data extension
Definition: muse_image.h:46
cpl_size muse_pixtable_get_nrow(const muse_pixtable *aPixtable)
get the number of rows within the pixel table
#define MUSE_HDR_PT_LHI
FITS header keyword contains the upper limit of the data in spectral direction.
#define MUSE_HDR_PT_ILO
FITS header keyword contains the lowest IFU number in the data.
int muse_pixtable_get_type(muse_pixtable *aPixtable)
Determine the type of pixel table.
#define MUSE_HDR_PT_EXP_LST
FITS header keyword defining the last row index for a given exposure.
Definition: muse_pixtable.h:98
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.
cpl_error_code muse_pixtable_erase_ifu_slice(muse_pixtable *aPixtable, unsigned char aIFU, unsigned short aSlice)
Erase pixel table rows related to one slice of one IFU.
cpl_error_code muse_pixtable_dump(muse_pixtable *aPixtable, cpl_size aStart, cpl_size aCount, unsigned char aDisplayHeader)
Dump a MUSE pixel table to the screen, resolving the origin column.
cpl_image * stat
the statistics extension
Definition: muse_image.h:64
#define MUSE_HDR_PT_FFCORR
#define MUSE_HDR_PT_TYPE
Pixel table "type" stored in the FITS header.
Definition: muse_pixtable.h:83
const char * muse_pfits_get_bunit(const cpl_propertylist *aHeaders)
find out the unit string
Definition: muse_pfits.c:223
#define MUSE_HDR_PT_MERGED
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_error_code muse_pixtable_origin_copy_offsets(muse_pixtable *aOut, muse_pixtable *aFrom, unsigned int aNum)
Copy MUSE_HDR_PT_IFU_SLICE_OFFSET keywords between pixel tables.
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.
#define MUSE_PIXTABLE_DATA
Definition: muse_pixtable.h:48
unsigned int muse_imagelist_get_size(muse_imagelist *aList)
Return the number of stored images.
muse_imagelist * muse_pixtable_to_imagelist(muse_pixtable *aPixtable)
Project a pixel table with data from one IFU back onto its image.
void muse_trace_polys_delete(cpl_polynomial *aPolys[])
Delete the multi-polynomial array created in relation to tracing.
#define MUSE_HDR_PT_SHI
FITS header keyword contains the highest slice number in the data.
cpl_image * dq
the data quality extension
Definition: muse_image.h:56
cpl_error_code muse_pixtable_append_ff(muse_pixtable *aPixtable, muse_image *aFF, cpl_table *aTrace, cpl_table *aWave, float aSampling)
Create flat-field spectrum and append to pixel table.
cpl_error_code muse_pixtable_restrict_wavelength(muse_pixtable *aPixtable, double aLow, double aHigh)
Restrict a pixel table to a certain wavelength range.
cpl_table * muse_cpltable_new(const muse_cpltable_def *aDef, cpl_size aLength)
Create an empty table according to the specified definition.
double muse_pfits_get_crpix(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS reference point
Definition: muse_pfits.c:401
cpl_boolean muse_pixtable_is_skysub(muse_pixtable *aPixtable)
Determine whether the pixel table is sky subtracted.
unsigned int muse_pixtable_origin_get_x(uint32_t aOrigin, muse_pixtable *aPixtable, cpl_size aRow)
Get the horizontal coordinate from the encoded 32bit origin number.
#define MUSE_PIXTABLE_WEIGHT
Definition: muse_pixtable.h:56
Structure definition of MUSE pixel table.
#define MUSE_HDR_PT_IFU_SLICE_OFFSET
FITS header keyword for the horizontal slice offset on the CCD.
Definition: muse_pixtable.h:74
#define MUSE_WCS_KEYS
regular expression for WCS properties
Definition: muse_wcs.h:48
muse_image * muse_imagelist_get(muse_imagelist *aList, unsigned int aIdx)
Get the muse_image of given list index.
cpl_array * muse_cpltable_extract_column(cpl_table *aTable, const char *aColumn)
Create an array from a section of a column.
#define MUSE_PIXTABLE_FF_EXT
Definition: muse_pixtable.h:59
muse_pixtable_wcs muse_pixtable_wcs_check(muse_pixtable *aPixtable)
Check the state of the world coordinate system of a pixel table.
#define MUSE_HDR_PT_YLO
FITS header keyword contains the lower limit of the data in y-direction.
#define MUSE_HDR_PT_SKYSUB
cpl_boolean muse_pixtable_is_fluxcal(muse_pixtable *aPixtable)
Determine whether the pixel table is flux calibrated.
#define MUSE_HDR_PT_ILLUM_REGEXP
cpl_table * ffspec
A flat-field spectrum.
cpl_array * muse_cplarray_interpolate_linear(const cpl_array *aTargetAbscissa, const cpl_array *aSourceAbscissa, const cpl_array *aSourceOrdinate)
Linear interpolation of a 1d array.
#define MUSE_HDR_PT_IHI
FITS header keyword contains the highest IFU number in the data.
unsigned int muse_pixtable_origin_get_offset(muse_pixtable *aPixtable, unsigned int aExpNum, unsigned short aIFU, unsigned short aSlice)
Get the slice offset from the pixel table header.
uint32_t muse_pixtable_origin_encode(unsigned int aX, unsigned int aY, unsigned short aIFU, unsigned short aSlice, unsigned int aOffset)
Encode the three CCD coordinates defining the origin of one MUSE pixel into a 32bit integer...
muse_pixtable * muse_pixtable_load_window(const char *aFilename, cpl_size aStart, cpl_size aNRows)
Load a range of rows from the table and all the FITS headers of a MUSE pixel table from a file...
muse_pixtable_operation
Type of operation to apply to a MUSE pixel table.
muse_pixtable * muse_pixtable_create(muse_image *aImage, cpl_table *aTrace, cpl_table *aWave, cpl_table *aGeoTable)
Create the pixel table for one CCD.
muse_pixtable * muse_pixtable_load_restricted_wavelength(const char *aFilename, double aLambdaMin, double aLambdaMax)
Load a pixel table from file and cut down the wavelength range.
#define MUSE_PIXTABLE_ORIGIN
Definition: muse_pixtable.h:54
cpl_error_code muse_cplpropertylist_update_long_long(cpl_propertylist *aHeader, const char *aKeyword, cpl_size aValue)
Update an integer-like property irrespective of the real type.
#define MUSE_HDR_PT_LLO
FITS header keyword contains the lower limit of the data in spectral direction.
cpl_error_code muse_pixtable_save(muse_pixtable *aPixtable, const char *aFilename)
Save a MUSE pixel table to a file on disk.
cpl_error_code muse_cpltable_copy_array(cpl_table *aTable, const char *aColumn, const cpl_array *aArray)
Copy an array into a table.
unsigned short muse_pixtable_origin_get_ifu(uint32_t aOrigin)
Get the IFU number from the encoded 32bit origin number.
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_pixtable_restrict_xpos(muse_pixtable *aPixtable, double aLo, double aHi)
Restrict a pixel table to a certain x coordinate range.
#define MUSE_PIXTABLE_STAT
Definition: muse_pixtable.h:50
cpl_error_code muse_pixtable_from_imagelist(muse_pixtable *aPixtable, muse_imagelist *aList)
Get pixel table values back from a per-IFU imagelist.
Handling of "mask" files.
Definition: muse_mask.h:43
muse_imagelist * muse_imagelist_new(void)
Create a new (empty) MUSE image list.
const muse_cpltable_def muse_pixtable_def[]
MUSE pixel table definition.
#define MUSE_PIXTABLE_LAMBDA
Definition: muse_pixtable.h:53
cpl_table * muse_geo_table_extract_ifu(const cpl_table *aTable, const unsigned char aIFU)
Extract the part of a geometry table dealing with a given IFU.
Definition: muse_geo.c:235
cpl_error_code muse_quadrants_coords_to_raw(cpl_propertylist *aHeader, int *aX, int *aY)
Convert coordinates of a trimmed image to raw-image coordinates.
#define MUSE_HDR_PT_SLO
FITS header keyword contains the lowest slice number in the data.
cpl_propertylist * header
the FITS header
Definition: muse_mask.h:56
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.
#define MUSE_HDR_PT_RVCORR
Definition of a cpl table structure.
muse_image * muse_image_new(void)
Allocate memory for a new muse_image object.
Definition: muse_image.c:66
cpl_table * muse_resampling_spectrum(muse_pixtable *aPixtable, double aBinwidth)
Resample the selected pixels of a pixel table into a spectrum.
unsigned int muse_pixtable_origin_get_y(uint32_t aOrigin)
Get the vertical coordinate from the encoded 32bit origin number.
cpl_mask * mask
The mask data.
Definition: muse_mask.h:49
#define MUSE_PIXTABLE_YPOS
Definition: muse_pixtable.h:52
void muse_pixtable_delete(muse_pixtable *aPixtable)
Deallocate memory associated to a pixel table object.
cpl_error_code muse_pixtable_restrict_ypos(muse_pixtable *aPixtable, double aLo, double aHi)
Restrict a pixel table to a certain y coordinate range.
#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.
static cpl_error_code muse_pixtable_fix_exp_headers(muse_pixtable *aPixtable)
Fix the exposure ranges in the header of a pixel table.
cpl_propertylist * header
The FITS header.
cpl_error_code muse_pixtable_flux_multiply(muse_pixtable *aPixtable, double aScale)
Scale the flux of a pixel table with correct treatment of variance.
#define MUSE_HDR_PT_XHI
FITS header keyword contains the upper limit of the data in x-ion.
cpl_boolean muse_pixtable_is_rvcorr(muse_pixtable *aPixtable)
Determine whether the pixel table is radial-velocity corrected.