X-shooter Pipeline Reference Manual 3.8.15
xsh_merge_ord.c
Go to the documentation of this file.
1/* *
2 * This file is part of the ESO X-shooter Pipeline *
3 * Copyright (C) 2006 European Southern Observatory *
4 * *
5 * This library is free software; you can redistribute it and/or modify *
6 * it under the terms of the GNU General Public License as published by *
7 * the Free Software Foundation; either version 2 of the License, or *
8 * (at your option) any later version. *
9 * *
10 * This program is distributed in the hope that it will be useful, *
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
13 * GNU General Public License for more details. *
14 * *
15 * You should have received a copy of the GNU General Public License *
16 * along with this program; if not, write to the Free Software *
17 * Foundation, 51 Franklin St, Fifth Floor, Boston, MA 02111-1307 USA *
18 * */
19
20/*
21 * $Author: amodigli $
22 * $Date: 2012-12-18 14:15:45 $
23 * $Revision: 1.71 $
24 */
25
26#ifdef HAVE_CONFIG_H
27#include <config.h>
28#endif
29
30/*----------------------------------------------------------------------------*/
37/*----------------------------------------------------------------------------*/
40/*-----------------------------------------------------------------------------
41 Includes
42 -----------------------------------------------------------------------------*/
43
44#include <math.h>
45#include <xsh_drl.h>
46
47#include <xsh_utils_table.h>
48#include <xsh_badpixelmap.h>
49#include <xsh_data_pre.h>
50#include <xsh_dfs.h>
51#include <xsh_pfits.h>
52#include <xsh_error.h>
53#include <xsh_msg.h>
54#include <xsh_fit.h>
55#include <xsh_data_instrument.h>
56#include <xsh_data_rec.h>
57#include <xsh_data_spectrum.h>
58#include <xsh_ifu_defs.h>
59#include <cpl.h>
60
61/*-----------------------------------------------------------------------------
62 Functions prototypes
63 ----------------------------------------------------------------------------*/
64/*-----------------------------------------------------------------------------
65 Implementation
66 ----------------------------------------------------------------------------*/
67
68/*****************************************************************************/
92/*****************************************************************************/
93static void xsh_merge_point( double flux_a, double weight_a,
94 double flux_b, double weight_b, double *flux_res,
95 double *err_res)
96{
97 XSH_ASSURE_NOT_NULL( flux_res);
98 XSH_ASSURE_NOT_NULL( err_res);
99
100 double tmp_val = 1.0/(weight_a+weight_b);
101 *flux_res = (weight_a*flux_a+weight_b*flux_b) * tmp_val;
102 *err_res = sqrt(tmp_val);
103
104 cleanup:
105 return;
106}
107#if 0
108static void
109xsh_spectrum_clean(xsh_spectrum *spectrum,const int decode_bp){
110
111 //double* spectrum_flux = xsh_spectrum_get_flux( spectrum);
112 //double* spectrum_errs = xsh_spectrum_get_errs( spectrum);
113 //int* spectrum_qual = xsh_spectrum_get_qual( spectrum);
114
115 cpl_mask* mask = xsh_qual_to_cpl_mask(spectrum->qual,decode_bp);
116
117 cpl_mask * filtered = xsh_bpm_filter(mask,5,5,CPL_FILTER_CLOSING);
118 xsh_free_mask(&mask);
119 cpl_image_save(spectrum->flux, "before.fits", CPL_TYPE_FLOAT, NULL, CPL_IO_CREATE);
120 cpl_image_reject_from_mask(spectrum->flux,filtered);
121
122 cpl_image_reject_from_mask(spectrum->errs,filtered);
123 cpl_detector_interpolate_rejected(spectrum->flux);
124 cpl_detector_interpolate_rejected(spectrum->errs);
125 cpl_image_save(spectrum->flux, "after.fits", CPL_TYPE_FLOAT, NULL, CPL_IO_CREATE);
126 xsh_free_mask(&filtered);
127
128 return;
129
130}
131#endif
132/*****************************************************************************/
147/*****************************************************************************/
148static cpl_frame * xsh_merge_ord_with_tag(cpl_frame *rec_frame,
149 xsh_instrument *instrument, int merge_par, const char *tag) {
150 cpl_frame *res_frame = NULL;
151 xsh_rec_list *rec_list = NULL;
152 xsh_spectrum *spectrum = NULL;
153 //int i;
154 int iorder = 0, islit = 0;
155 int spectrum_size_lambda = 0;
156 int spectrum_size_slit = 0;
157 double lambda_min = 0.0, lambda_max = 0.0, lambda_step = 0.0;
158 double *spectrum_flux = NULL;
159 double *spectrum_errs = NULL;
160 int *spectrum_qual = NULL;
161 int spectrum_flux_index = 0;
162 int *spectrum_by_lambda = NULL;
163 const char* rec_pcatg = NULL;
164 char* spectrum_name = NULL;
165 const char* name = NULL;
166 cpl_propertylist* header = NULL;
167 int naxis = 0;
168 int eso_mode = 0;
169 int decode_bp = instrument->decode_bp;
170 /* Check parameters */
171 XSH_ASSURE_NOT_NULL( rec_frame);
174
175
176 /* Load rectified table */
177 name = cpl_frame_get_filename(rec_frame);
178 header = cpl_propertylist_load(name, 0);
179 naxis = xsh_pfits_get_naxis(header);
180 xsh_free_propertylist(&header);
181 if (naxis == 2) {
182 check( rec_list = xsh_rec_list_load_eso( rec_frame, instrument));
183 eso_mode = 1;
184 } else if (naxis == 1) {
185 check( rec_list = xsh_rec_list_load_eso_1d( rec_frame, instrument));
186 eso_mode = 1;
187 } else {
188 check( rec_list = xsh_rec_list_load( rec_frame, instrument));
189 }
190
191 check(lambda_min = xsh_pfits_get_rectify_lambda_min( rec_list->header));
192 check( lambda_max = xsh_pfits_get_rectify_lambda_max( rec_list->header));
194
195 /* Product spectrum */
196 check( rec_pcatg = xsh_pfits_get_pcatg( rec_list->header));
197
198 /* Check if 1D or 2D spectrum */
199 if (strstr(rec_pcatg, XSH_ORDER2D) != NULL
200 || strstr(rec_pcatg, XSH_MERGE2D) != NULL
201 || strstr(rec_pcatg, XSH_COMBINED_OFFSET_2D_SLIT) != NULL
202 || strstr(rec_pcatg, XSH_COMBINED_OFFSET_2D_IFU) != NULL) {
203
204 double slit_min = 0.0, slit_max = 0.0, slit_step = 0.0;
205
206 check( slit_min = xsh_pfits_get_rectify_space_min( rec_list->header));
207 check( slit_max = xsh_pfits_get_rectify_space_max( rec_list->header));
209
211 "Construct 2D spectrum %f %f %f", slit_min, slit_max, slit_step);
213 "lambda min %f max %f step %f", lambda_min, lambda_max, lambda_step);
214 check(
215 spectrum = xsh_spectrum_2D_create( lambda_min, lambda_max, lambda_step, slit_min, slit_max, slit_step));
216 } else {
217 xsh_msg_dbg_high("Construct 1D spectrum");
219 "lambda min %f max %f step %f", lambda_min, lambda_max, lambda_step);
220 check(
221 spectrum = xsh_spectrum_1D_create( lambda_min, lambda_max, lambda_step));
222 }
223
224 spectrum_name = xsh_stringcat_any(tag, ".fits", (void*)NULL);
225
226 check( spectrum_size_lambda = xsh_spectrum_get_size_lambda( spectrum));
227 check( spectrum_size_slit = xsh_spectrum_get_size_slit( spectrum));
228 check( spectrum_flux = xsh_spectrum_get_flux( spectrum));
229 check( spectrum_errs = xsh_spectrum_get_errs( spectrum));
230 check( spectrum_qual = xsh_spectrum_get_qual( spectrum));
231 XSH_CALLOC( spectrum_by_lambda, int, spectrum->size);
232 //xsh_msg("spectrum_size_slit=%d",spectrum_size_slit);
233 xsh_msg_dbg_medium ( "Spectrum size (Lambda X slit) = (%d X %d) \n"
234 "Lambda %f to %f", spectrum_size_lambda, spectrum_size_slit,
235 lambda_min, lambda_min+lambda_step*(spectrum_size_lambda-1));
236 /*
237 xsh_msg("ord_max=%d spectrum_size_slit %d",
238 rec_list->size-1,spectrum_size_slit);
239 */
240 int problem_slit_size = 0;
241 /* loop over spatial direction: we merge each spectrum along the slit */
242 for (islit = 0; islit < spectrum_size_slit; islit++) {
243 for (iorder = rec_list->size - 1; iorder >= 0; iorder--) {
244 int nlambda = xsh_rec_list_get_nlambda(rec_list, iorder);
245 double* lambda = xsh_rec_list_get_lambda(rec_list, iorder);
246 float *flux = xsh_rec_list_get_data1(rec_list, iorder);
247 float *errs = xsh_rec_list_get_errs1(rec_list, iorder);
248 int *qual = xsh_rec_list_get_qual1(rec_list, iorder);
249 int nslit = xsh_rec_list_get_nslit(rec_list, iorder);
250 int ilambda = 0;
251 double n;
252 /*
253 xsh_msg("iorder=%d islit %d nslit=%d nlambda=%d",
254 iorder,islit,nslit,nlambda);
255 */
256 if (spectrum_size_slit != nslit) {
257 /* AMO: REFLEX: move this check out of the loop to decrease log size
258 xsh_msg_warning("spectrum size slit %d is different from order size %d",
259 spectrum_size_slit, nslit);
260 */
261 problem_slit_size = 1;
262 }
263
265 "slit %d lambda %f %f", islit, lambda[0], lambda[nlambda-1]);
266 /*
267 xsh_msg( "slit %d lambda %f %f", islit,
268 lambda[0], lambda[nlambda-1]);
269 */
270 n = (lambda[0] - lambda_min) / lambda_step;
271 spectrum_flux_index = (int) xsh_round_double(n);
272
274 "iorder %d begin at index %f --> %d", iorder, n, spectrum_flux_index);
275
276 if (islit < nslit) {
277 int ilambda_off=spectrum_flux_index+ islit*spectrum_size_lambda;
278 int islit_nlambda = islit * nlambda;
279 /* look over wavelength */
280 for (ilambda = 0; ilambda < nlambda; ilambda++) {
281
282
283 int i_rec_order = ilambda + islit_nlambda;
284 int i_spectrum = ilambda_off + ilambda;
285 /* skip merging if error associated to given point is zero.
286 * this may happen if extraction has not found enough slices
287 * with good points to get the cross-order spectrum extraction profile
288 * In this case the merged flux and associated error os 0, but the
289 * qualifier is the same as the one of the original spectrum
290 */
291 if (errs[i_rec_order] == 0) {
292 spectrum_qual[i_spectrum] = qual[i_rec_order];
293 continue;
294 }
295
296 /* The following setting may be useful for debug purpose to check
297 * at which wavelength we are
298 * double wav=ilambda*lambda_step+lambda[0];
299 */
300 /* initialise the merged spectrum value with the corresponding values
301 * from the extracted profile. Uses the array spectrum_by_lambda to
302 * check if such event already occurred for a If this is the case at
303 * a given wavelength, it means that there is overlap between the
304 * current order spectrum data point and the one of the previous
305 * order. In that case we need to merge the intensity of the two data
306 * points. Correspondingly we distinguish three cases:
307 * A good, B bad: use only A contribute
308 * A bad, B good: use only B contribute
309 * A, B good: use weighted mean from A and B contributes
310 * Finally we set the result to the proper value.
311 */
312 if (spectrum_by_lambda[i_spectrum] == 0) {
313 spectrum_by_lambda[i_spectrum] = 1;
314 spectrum_flux[i_spectrum] = flux[i_rec_order];
315 spectrum_errs[i_spectrum] = errs[i_rec_order];
316 spectrum_qual[i_spectrum] = qual[i_rec_order];
317 } else {
318 double flux_a, err_a, weight_a;
319 double flux_b, err_b, weight_b;
320 double flux_res = 0.0, err_res = 1.0;
321 int qual_a, qual_b, qual_res;
322
323 flux_a = spectrum_flux[i_spectrum];
324 err_a = spectrum_errs[i_spectrum];
325 qual_a = spectrum_qual[i_spectrum];
326
327 flux_b = flux[i_rec_order];
328 err_b = errs[i_rec_order];
329 qual_b = qual[i_rec_order];
330
331 /* Predefine outputs as from first A spectra */
332 flux_res = flux_a;
333 err_res = err_a;
334 qual_res = qual_a;
335 //xsh_msg("wave=%g",lambda[ilambda]);
336 /* A good, B bad */
337 if (((qual_a & decode_bp) == 0) && ((qual_b & decode_bp) > 0)) {
338 flux_res = flux_a;
339 err_res = err_a;
340 qual_res = qual_a;
341
342 /* A bad, B good */
343 } else if (((qual_a & decode_bp) > 0) && ((qual_b & decode_bp) == 0)) {
344 flux_res = flux_b;
345 err_res = err_b;
346 qual_res = qual_b;
347
348 /* A and B good, or A and B bad */
349 } else {
350
351
352 weight_a = 1.0 / (err_a * err_a);
353 weight_b = 1.0 / (err_b * err_b);
354
355 check(xsh_merge_point( flux_a, weight_a, flux_b, weight_b, &flux_res, &err_res));
356 qual_res = qual_a | qual_b;
357
358 }
359
360 /* Save values */
361 spectrum_flux[i_spectrum] = flux_res;
362 spectrum_errs[i_spectrum] = err_res;
363 spectrum_qual[i_spectrum] = qual_res;
364 spectrum_by_lambda[i_spectrum] = 1;
365 }
366 //xsh_msg("loop lambda ok");
367 }
368 //xsh_msg("check slit ok");
369 }
370 }/* end order*/
371
372 }/*end slit*/
373 if (problem_slit_size == 1) {
374 xsh_msg_warning("spectrum size slit is different from order size");
375 }
376 /* Set the header from input frame */
377 if (eso_mode) {
378 cpl_propertylist_erase_regexp(rec_list->header, "^CRVAL1", 0);
379 cpl_propertylist_erase_regexp(rec_list->header, "^CRVAL2", 0);
380 cpl_propertylist_erase_regexp(rec_list->header, "^CRPIX1", 0);
381 cpl_propertylist_erase_regexp(rec_list->header, "^CRPIX2", 0);
382 cpl_propertylist_erase_regexp(rec_list->header, "^CDELT1", 0);
383 cpl_propertylist_erase_regexp(rec_list->header, "^CDELT2", 0);
384 }
385 check( cpl_propertylist_append( spectrum->flux_header, rec_list->header));
386
387 const char* bunit = xsh_pfits_get_bunit(spectrum->flux_header);
388 xsh_pfits_set_bunit(spectrum->errs_header, bunit);
389 /* AMO to verify quality of extraction */
390 //xsh_spectrum_clean(spectrum,decode_bp);
391
392 check( res_frame = xsh_spectrum_save(spectrum, spectrum_name,tag));
393 check( cpl_frame_set_tag( res_frame, tag));
394
395 cleanup: xsh_spectrum_free(&spectrum);
396 xsh_rec_list_free(&rec_list);
397 XSH_FREE( spectrum_by_lambda);
398 XSH_FREE( spectrum_name);
399 return res_frame;
400}
401/*****************************************************************************/
402/*****************************************************************************/
417/*****************************************************************************/
418cpl_frame * xsh_merge_ord( cpl_frame * rec_frame, xsh_instrument* instrument,
419 int merge_par,const char* rec_prefix)
420{
421 cpl_frame *res_frame = NULL;
422 xsh_msg("Merge slit orders");
423 check( res_frame = xsh_merge_ord_slitlet( rec_frame, instrument, merge_par,
424 CENTER_SLIT,rec_prefix));
425
426 cleanup:
427 return res_frame;
428}
429/*****************************************************************************/
446/*****************************************************************************/
447cpl_frame * xsh_merge_ord_slitlet( cpl_frame * rec_frame,
449 int merge_par, int slitlet,
450 const char* rec_prefix)
451{
452 cpl_frame *res_frame = NULL;
453 cpl_propertylist *header = NULL;
454 const char *rec_pcatg = NULL;
455 const char *tag_suf = NULL;
456 const char *filename = NULL;
457 char tag[256];
458
459 /*Check parameters */
460 XSH_ASSURE_NOT_NULL( rec_frame);
462
463 /* Check 1D or 2D spectrum */
464 check( filename = cpl_frame_get_filename( rec_frame));
465 check( header = cpl_propertylist_load( filename, 0));
466 check( rec_pcatg = xsh_pfits_get_pcatg( header));
467 check(tag_suf = cpl_frame_get_tag(rec_frame));
468
469 if ( strstr( rec_pcatg, XSH_FLUX_ORDER2D ) != NULL ) {
471 }
472 else if ( strstr( rec_pcatg, XSH_NORM_ORDER2D ) != NULL ) {
474 }
475// else if ( strstr( rec_pcatg, XSH_SKY_ORDER2D ) != NULL ) {
476// tag_suf = XSH_GET_TAG_FROM_SLITLET( XSH_SKY_MERGE2D, slitlet,
477// instrument);
478// }
479 else if ( strstr( rec_pcatg, XSH_ORDER2D ) != NULL ) {
480 tag_suf = XSH_GET_TAG_FROM_SLITLET( XSH_MERGE2D, slitlet,
481 instrument);
482 }
483 else if ( strstr( rec_pcatg, XSH_COMBINED_OFFSET_2D_SLIT ) != NULL){
484 /* No IFU */
486 instrument);
487 }
488 else if ( strstr( rec_pcatg, XSH_ORDER_EXT1D ) != NULL){
489 /* No IFU */
491 instrument);
492 }
493 else if ( strstr( rec_pcatg, XSH_ORDER_OXT1D ) != NULL){
494 /* No IFU */
496 instrument);
497 }
498 else if ( strstr( rec_pcatg, XSH_FLUX_OXT1D ) != NULL){
500 }
501 else if ( strstr( rec_pcatg, XSH_FLUX_ORDER1D ) != NULL){
503 }
504 else if ( strstr( rec_pcatg, XSH_NORM_ORDER1D ) != NULL){
506 }
507 else if ( strstr( rec_pcatg, XSH_ORDER1D ) != NULL){
508 tag_suf = XSH_GET_TAG_FROM_SLITLET( XSH_MERGE1D, slitlet,
509 instrument);
510 }
511 else{
512 tag_suf="";
513 xsh_msg_error("Invalid PRO.CATG %s for file %s", rec_pcatg, filename);
514 }
515 sprintf(tag,"%s_%s",rec_prefix,tag_suf);
516 check( res_frame = xsh_merge_ord_with_tag( rec_frame, instrument,
517 merge_par, tag));
518
519 cleanup:
520 xsh_free_propertylist( &header);
521 return res_frame ;
522}
523
524/*****************************************************************************/
525static void xsh_frame_set_shift_ref( cpl_frame *rec_frame, cpl_frame *shift_frame)
526{
527 cpl_propertylist *rec_header = NULL;
528 cpl_propertylist *shift_header = NULL;
529 const char *rec_name = NULL;
530 const char *shift_name = NULL;
531 double waveref, slitref;
532
533 XSH_ASSURE_NOT_NULL( rec_frame);
534 XSH_ASSURE_NOT_NULL( shift_frame);
535
536 check( rec_name = cpl_frame_get_filename( rec_frame));
537 check( shift_name = cpl_frame_get_filename( shift_frame));
538
539 check( rec_header = cpl_propertylist_load( rec_name, 0));
540 check( shift_header = cpl_propertylist_load( shift_name, 0));
541
542 waveref = xsh_pfits_get_shiftifu_lambdaref( shift_header);
543 slitref = xsh_pfits_get_shiftifu_slitref( shift_header);
544
545 if (cpl_error_get_code() == CPL_ERROR_NONE){
546 check( xsh_pfits_set_shiftifu_lambdaref( rec_header, waveref));
547 check( xsh_pfits_set_shiftifu_slitref( rec_header, slitref));
548 check( cpl_propertylist_save( rec_header, rec_name, CPL_IO_EXTEND));
549 }
551
552 cleanup:
553 xsh_free_propertylist( &rec_header);
554 xsh_free_propertylist( &shift_header);
555 return;
556}
557/*****************************************************************************/
558
559/*****************************************************************************/
574/*****************************************************************************/
575cpl_frameset* xsh_merge_ord_ifu(cpl_frameset* rec_frameset,
577 int merge_par,
578 const char* rec_prefix)
579{
580 cpl_frameset *result_set = NULL;
581 cpl_frameset *drl_frameset = NULL;
582 int slitlet;
583 int i=0;
584
585 XSH_ASSURE_NOT_NULL( rec_frameset);
587
588 xsh_msg("Merge IFU orders");
589 check( result_set = cpl_frameset_new());
590 check( drl_frameset = xsh_frameset_drl_frames( rec_frameset));
591
592 /* Loop over the 3 IFU slitlets */
593 for( slitlet = LOWER_IFU_SLITLET; slitlet <= UPPER_IFU_SLITLET; slitlet++){
594 cpl_frame * rec_frame = NULL ;
595 cpl_frame * ord_frame = NULL ;
596
597 check( rec_frame = cpl_frameset_get_frame( drl_frameset, i));
598
599 /* Now merge order for this slitlet */
600 check( ord_frame = xsh_merge_ord_slitlet( rec_frame, instrument,
601 merge_par,slitlet,rec_prefix));
602
603 check( xsh_frame_set_shift_ref( rec_frame, ord_frame));
604 xsh_msg( "Merge for Slitlet %s, %s", SlitletName[slitlet],
605 cpl_frame_get_filename( ord_frame)) ;
606 check( cpl_frameset_insert( result_set, ord_frame));
607 i++;
608 }
609
610 cleanup:
611 if (cpl_error_get_code() != CPL_ERROR_NONE){
612 xsh_free_frameset( &result_set);
613 }
614 xsh_free_frameset( &drl_frameset);
615 return result_set ;
616}
617/*--------------------------------------------------------------------------*/
static xsh_instrument * instrument
static float slit_step
static float lambda_step
cpl_mask * xsh_qual_to_cpl_mask(cpl_image *qual, const int decode_bp)
cpl_mask * xsh_bpm_filter(const cpl_mask *input_mask, cpl_size kernel_nx, cpl_size kernel_ny, cpl_filter_mode filter)
Allows the growing and shrinking of bad pixel masks. It can be used to e.g. set pixels to bad if the ...
int xsh_rec_list_get_nslit(xsh_rec_list *list, int idx)
Definition: xsh_data_rec.c:827
double * xsh_rec_list_get_lambda(xsh_rec_list *list, int idx)
xsh_rec_list * xsh_rec_list_load(cpl_frame *frame, xsh_instrument *instrument)
load an rec list from a frame
Definition: xsh_data_rec.c:289
float * xsh_rec_list_get_data1(xsh_rec_list *list, int idx)
float * xsh_rec_list_get_errs1(xsh_rec_list *list, int idx)
xsh_rec_list * xsh_rec_list_load_eso(cpl_frame *frame, xsh_instrument *instrument)
Definition: xsh_data_rec.c:407
xsh_rec_list * xsh_rec_list_load_eso_1d(cpl_frame *frame, xsh_instrument *instrument)
Definition: xsh_data_rec.c:563
int * xsh_rec_list_get_qual1(xsh_rec_list *list, int idx)
int xsh_rec_list_get_nlambda(xsh_rec_list *list, int idx)
Definition: xsh_data_rec.c:865
void xsh_rec_list_free(xsh_rec_list **list)
free memory associated to a rec_list
Definition: xsh_data_rec.c:760
int xsh_spectrum_get_size_lambda(xsh_spectrum *s)
Get lambda axis size of spectrum.
double * xsh_spectrum_get_errs(xsh_spectrum *s)
Get errs of spectrum.
int * xsh_spectrum_get_qual(xsh_spectrum *s)
Get qual of spectrum.
cpl_frame * xsh_spectrum_save(xsh_spectrum *s, const char *filename, const char *tag)
save a spectrum
double * xsh_spectrum_get_flux(xsh_spectrum *s)
Get flux of spectrum.
xsh_spectrum * xsh_spectrum_2D_create(double lambda_min, double lambda_max, double lambda_step, double slit_min, double slit_max, double slit_step)
Create a 2D spectrum structure.
xsh_spectrum * xsh_spectrum_1D_create(double lambda_min, double lambda_max, double lambda_step)
Create a 1D spectrum structure.
void xsh_spectrum_free(xsh_spectrum **s)
free memory associated to an 1D spectrum
int xsh_spectrum_get_size_slit(xsh_spectrum *s)
Get slit axis ize of spectrum.
#define check(COMMAND)
Definition: xsh_error.h:71
#define xsh_error_reset()
Definition: xsh_error.h:87
#define XSH_ASSURE_NOT_NULL(pointer)
Definition: xsh_error.h:99
cpl_frameset * xsh_merge_ord_ifu(cpl_frameset *rec_frameset, xsh_instrument *instrument, int merge_par, const char *rec_prefix)
Merge orders of the rectified frame using merge parameters.
cpl_frame * xsh_merge_ord(cpl_frame *rec_frame, xsh_instrument *instrument, int merge_par, const char *rec_prefix)
Merge orders of the rectified frame using merge parameters.
static void xsh_merge_point(double flux_a, double weight_a, double flux_b, double weight_b, double *flux_res, double *err_res)
Compute flux and error associated to each merged spectrum point.
Definition: xsh_merge_ord.c:93
static void xsh_frame_set_shift_ref(cpl_frame *rec_frame, cpl_frame *shift_frame)
cpl_frame * xsh_merge_ord_slitlet(cpl_frame *rec_frame, xsh_instrument *instrument, int merge_par, int slitlet, const char *rec_prefix)
Merge orders of the rectified frame using merge parameters.
static cpl_frame * xsh_merge_ord_with_tag(cpl_frame *rec_frame, xsh_instrument *instrument, int merge_par, const char *tag)
Merge the orders.
#define xsh_msg_warning(...)
Print an warning message.
Definition: xsh_msg.h:88
#define xsh_msg_dbg_medium(...)
Definition: xsh_msg.h:44
#define xsh_msg_error(...)
Print an error message.
Definition: xsh_msg.h:62
#define xsh_msg(...)
Print a message on info level.
Definition: xsh_msg.h:121
#define xsh_msg_dbg_high(...)
Definition: xsh_msg.h:40
double xsh_pfits_get_rectify_space_max(cpl_propertylist *plist)
find out the rectify SPACE max
Definition: xsh_pfits.c:3370
const char * xsh_pfits_get_bunit(const cpl_propertylist *plist)
find out the BUNIT
Definition: xsh_pfits.c:1467
double xsh_pfits_get_rectify_lambda_min(cpl_propertylist *plist)
find out the rectify lambda min
Definition: xsh_pfits.c:3319
double xsh_pfits_get_rectify_space_min(cpl_propertylist *plist)
find out the rectify space min
Definition: xsh_pfits.c:3353
double xsh_pfits_get_shiftifu_lambdaref(cpl_propertylist *plist)
Definition: xsh_pfits.c:4030
double xsh_pfits_get_rectify_bin_lambda(cpl_propertylist *plist)
find out the rectify lambda binning
Definition: xsh_pfits.c:3284
void xsh_pfits_set_shiftifu_lambdaref(cpl_propertylist *plist, double value)
Definition: xsh_pfits.c:4015
double xsh_pfits_get_rectify_bin_space(cpl_propertylist *plist)
find out the rectify space (slit) binning
Definition: xsh_pfits.c:3302
const char * xsh_pfits_get_pcatg(const cpl_propertylist *plist)
find out the pcatg
Definition: xsh_pfits.c:1627
void xsh_pfits_set_shiftifu_slitref(cpl_propertylist *plist, double value)
Definition: xsh_pfits.c:4044
double xsh_pfits_get_shiftifu_slitref(cpl_propertylist *plist)
Definition: xsh_pfits.c:4104
void xsh_pfits_set_bunit(cpl_propertylist *plist, const char *value)
Write the BUNIT value.
Definition: xsh_pfits.c:2606
int xsh_pfits_get_naxis(const cpl_propertylist *plist)
find out the NAXIS value
Definition: xsh_pfits.c:209
double xsh_pfits_get_rectify_lambda_max(cpl_propertylist *plist)
find out the rectify lambda max
Definition: xsh_pfits.c:3336
void xsh_free_frameset(cpl_frameset **f)
Deallocate a frame set and set the pointer to NULL.
Definition: xsh_utils.c:2254
void xsh_free_mask(cpl_mask **m)
Deallocate an image mask and set the pointer to NULL.
Definition: xsh_utils.c:2149
char * xsh_stringcat_any(const char *s,...)
Concatenate an arbitrary number of strings.
Definition: xsh_utils.c:1925
long xsh_round_double(double x)
Computes round(x)
Definition: xsh_utils.c:4383
void xsh_free_propertylist(cpl_propertylist **p)
Deallocate a property list and set the pointer to NULL.
Definition: xsh_utils.c:2179
cpl_propertylist * header
Definition: xsh_data_rec.h:89
cpl_propertylist * errs_header
cpl_propertylist * flux_header
cpl_image * flux
cpl_image * qual
cpl_image * errs
int n
Definition: xsh_detmon_lg.c:92
cpl_frameset * xsh_frameset_drl_frames(cpl_frameset *set)
extract DRL specific frames from frameset
Definition: xsh_dfs.c:319
#define XSH_NORM_ORDER2D
Definition: xsh_dfs.h:986
#define XSH_ORDER_EXT1D
Definition: xsh_dfs.h:1029
#define XSH_FLUX_MOXT1D
Definition: xsh_dfs.h:1052
#define XSH_FLUX_MERGE2D
Definition: xsh_dfs.h:1024
#define XSH_MERGE_OXT1D
Definition: xsh_dfs.h:1047
#define XSH_MERGE_EXT1D
Definition: xsh_dfs.h:593
#define XSH_GET_TAG_FROM_SLITLET(TAG, slitlet, instr)
Definition: xsh_dfs.h:1559
#define XSH_NORM_ORDER1D
Definition: xsh_dfs.h:981
#define XSH_ORDER2D
Definition: xsh_dfs.h:588
#define XSH_ORDER1D
Definition: xsh_dfs.h:587
#define XSH_NORM_MERGE1D
Definition: xsh_dfs.h:991
#define XSH_COMBINED_OFFSET_2D_IFU
Definition: xsh_dfs.h:838
#define XSH_FLUX_ORDER2D
Definition: xsh_dfs.h:1014
#define XSH_GET_TAG_FROM_ARM(TAG, instr)
Definition: xsh_dfs.h:1548
#define XSH_NORM_MERGE2D
Definition: xsh_dfs.h:997
#define XSH_COMBINED_OFFSET_2D_SLIT
Definition: xsh_dfs.h:854
#define XSH_MERGE2D
Definition: xsh_dfs.h:594
#define XSH_MERGE1D
Definition: xsh_dfs.h:592
#define XSH_ORDER_OXT1D
Definition: xsh_dfs.h:1034
#define XSH_FLUX_MERGE1D
Definition: xsh_dfs.h:1019
#define XSH_FLUX_ORDER1D
Definition: xsh_dfs.h:1003
#define XSH_FLUX_OXT1D
Definition: xsh_dfs.h:1009
@ LOWER_IFU_SLITLET
Definition: xsh_ifu_defs.h:42
@ UPPER_IFU_SLITLET
Definition: xsh_ifu_defs.h:42
@ CENTER_SLIT
Definition: xsh_ifu_defs.h:41
static const char * SlitletName[]
Definition: xsh_ifu_defs.h:48
#define XSH_FREE(POINTER)
Definition: xsh_utils.h:92
#define XSH_CALLOC(POINTER, TYPE, SIZE)
Definition: xsh_utils.h:56