X-shooter Pipeline Reference Manual 3.8.15
xsh_rectify.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 F1ITNESS 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.151 $
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_baryvel.h>
48#include <xsh_badpixelmap.h>
49#include <xsh_data_rec.h>
51#include <xsh_data_pre.h>
52#include <xsh_data_order.h>
53#include <xsh_data_wavesol.h>
55#include <xsh_dfs.h>
56#include <xsh_pfits.h>
57#include <xsh_error.h>
58#include <xsh_msg.h>
59#include <xsh_fit.h>
60#include <xsh_model_io.h>
61#include <xsh_model_kernel.h>
62#include <xsh_ifu_defs.h>
63#include <xsh_data_dispersol.h>
64#include <cpl.h>
65#include <xsh_utils.h>
66#include <xsh_utils_ifu.h>
67#include <xsh_utils_table.h>
68#include <xsh_rectify.h>
69#include <xsh_model_utils.h>
70
71#define new_bp 1
72/*#include <omp.h> */
73
74/*-----------------------------------------------------------------------------
75 Functions prototypes
76 -----------------------------------------------------------------------------*/
77
78#define SLIT_FRAC_STEP 50
79
80static cpl_vector * rec_profile = NULL ;
81static cpl_vector * err_profile = NULL ;
82
83static double compute_shift_with_localization( cpl_frame *loc_frame,
84 cpl_frame *loc0_frame);
85
86static double compute_shift_with_kw( cpl_propertylist * header,
87 xsh_rectify_param *rectify_par,
88 double **ref_ra, double **ref_dec, int flag);
89
90static cpl_frame*
91xsh_shift( cpl_frame *rec_frame, xsh_instrument *instrument,
92 const char *fname, double slit_shift, cpl_frame** res_frame_ext);
93
94static void xsh_frame_set_shiftifu_ref( cpl_propertylist *header,
95 cpl_frame *shift_frame);
96
97/*-----------------------------------------------------------------------------
98 Implementation
99 -----------------------------------------------------------------------------*/
100/*
101static void get_errors( float *err, int ix, int iy, int nx, int ny,
102 double radius, float *err_val, int* qual_val)
103{
104 int xmin, xmax, ymin, ymax ;
105 int i, j, x, y, n ;
106 double sum2 ;
107 double *err_vect = NULL;
108 int err_vect_size = 0;
109
110 XSH_ASSURE_NOT_NULL( err);
111 XSH_ASSURE_NOT_NULL( err_val);
112 XSH_ASSURE_NOT_NULL( qual_val);
113
114 // Get errors around ix, iy
115 // Sum squares
116 // Divide by N
117
118 n = 0 ;
119 sum2 = 0 ;
120 xmin = ix - radius ;
121 xmax = ix + radius ;
122 ymin = iy - radius ;
123 ymax = iy + radius ;
124
125
126 err_vect_size = (2*radius+1)*(2*radius+1);
127
128 XSH_CALLOC( err_vect, double, err_vect_size);
129
130 for( y = ymin ; y < ymax ; y++ ){
131 for( x = xmin ; x < xmax ; x++ ) {
132 if ( y < 0 || y >= ny || x < 0 || x >= nx ) continue;
133 err_vect[n] = err[x+y*nx];
134 n++;
135 }
136 }
137
138 for(i=0; i< n; i++){
139 for(j=i; j< n; j++){
140 if (i==j){
141 sum2 += err_vect[i]*err_vect[j];
142 }
143 else{
144 sum2 += 2*err_vect[i]*err_vect[j];
145 }
146 }
147 }
148
149 if (n != 0){
150 *err_val = sqrt( sum2 )/(double)n ;
151 }
152 else{
153 *err_val = 1;
154 *qual_val = QFLAG_OUTSIDE_DATA_RANGE;
155 }
156 cleanup:
157 XSH_FREE( err_vect);
158}
159
160*/
161/*****************************************************************************/
162#if 0
163static float
164xsh_fill_gauss(const int j, const float A, const float sigma, const float bkg,
165 const float y_c)
166{
167
168 /* add Gauss profile positioned where the trace is
169 * I=A*exp( -0.5*( (X-Xo)/sigma )^2 ) +bkg
170 */
171 float gauss;
172 float arg = 0;
173 float y;
174
175 /* uniform along X (horizontal direction) Gauss along Y (vertical) */
176
177 y=j;
178
179 arg = (y - y_c) / sigma;
180 arg *= arg;
181 gauss = (A * expf(-0.5 * arg) + bkg);
182 //xsh_msg("y=%g y_c=%g arg=%g gauss=%g",y,y_c,arg,gauss);
183 return gauss;
184}
185#endif
186/*****************************************************************************/
187static void
188xsh_rec_list_rectify( xsh_rec_list *rec_list, int iorder, int irec,
189 xsh_pre *sci_pre,cpl_image* var_img,
190 double fx, double fy, double radius, cpl_vector *profile,
191 cpl_vector* profile2,
192 double mult)
193{
194 int nx,ny;
195 float *rec_flux = NULL;
196 float *rec_errs = NULL;
197 int *rec_qual = NULL;
198 int ix,iy;
199
200 int xfirst=0;
201 int xlast=0;
202
203 int yfirst=0;
204 int ylast=0;
205
206 register int i=0;
207 register int j=0;
208 register int j_nx=0;
209
210 XSH_ASSURE_NOT_NULL( rec_list);
211 XSH_ASSURE_NOT_NULL( sci_pre);
212 XSH_ASSURE_NOT_NULL( profile);
213
214 check( rec_flux = xsh_rec_list_get_data1( rec_list, iorder));
215 check( rec_errs = xsh_rec_list_get_errs1( rec_list, iorder));
216 check( rec_qual = xsh_rec_list_get_qual1( rec_list, iorder));
217
218 check( nx = xsh_pre_get_nx( sci_pre));
219 check( ny = xsh_pre_get_ny( sci_pre));
220
221 ix = floor( fx ) ;
222 iy = floor( fy ) ;
223 xfirst=ceil(fx-radius);
224 xlast=floor(fx+radius);
225
226 yfirst=ceil(fy-radius);
227 ylast=floor(fy+radius);
228 /* to prevent to go outside the image we do the following check */
229 xfirst = (xfirst > 0) ? xfirst : 0;
230 xlast = (xlast < nx) ? xlast : nx;
231 yfirst = (yfirst > 0) ? yfirst : 0;
232 ylast = (ylast < ny) ? ylast : ny;
233
234
235 /* If outside detector */
236 if ( ix <0 || ix >= nx || iy < 0 || iy >= ny ) {
237 xsh_msg_dbg_high( "OUT_OF_RANGE at %d,%d", ix, iy);
238 rec_flux[irec] = 0;
239 rec_errs[irec] = 1;
240 rec_qual[irec] = QFLAG_OUTSIDE_DATA_RANGE;
241 }
242
243 /* Otherwise if inside detector */
244 else{
245 cpl_image *sci_img = NULL;
246 cpl_image *err_img = NULL;
247 cpl_image *qua_img = NULL;
248 float *data = NULL;
249 float *errs = NULL;
250 int *qual = NULL;
251
252 check( sci_img = xsh_pre_get_data( sci_pre));
253 check( err_img = xsh_pre_get_errs( sci_pre));
254 check( qua_img = xsh_pre_get_qual( sci_pre));
255
256 check( data = cpl_image_get_data_float( sci_img));
257 check( errs = cpl_image_get_data_float( err_img));
258 check( qual = cpl_image_get_data_int( qua_img));
259
260 int ix_iynx=ix+iy*nx;
261
262 /* If the kernel radius is zero */
263 if ( radius == 0.){
264 rec_flux[irec] = data[ix_iynx];
265 rec_errs[irec] = errs[ix_iynx];
266 rec_qual[irec] = qual[ix_iynx];
267 }
268
269 /* Otherwise if kernel radius is greater than zero */
270 else {
271
272#if new_bp
273 double flux_val=0, confidence=0 ;
274 check( flux_val = cpl_image_get_interpolated( sci_img, fx, fy,
275 profile, radius, profile, radius, &confidence));
276
277 /* Now propagate errors */
278 float err_val = 0;
279 check( err_val = cpl_image_get_interpolated( var_img, fx, fy,
280 profile2, radius, profile2, radius, &confidence));
281
282 /* Now propagate flags */
283 int qual_val = 0;
284 for(j=yfirst;j<ylast-1;j++) {
285 j_nx=j*nx;
286 for(i=xfirst;i<xlast-1;i++) {
287 qual_val |= qual[j_nx+i];
288 }
289 }
290
291 /* Save the interpolated values */
292 rec_flux[irec] = flux_val;
293 rec_errs[irec] = sqrt(err_val);
294 if (confidence < 1) {
295 /* If error or undefined or pixel within kernel radius of edge */
296 rec_qual[irec] = qual_val || QFLAG_MISSING_DATA;
297 } else {
298 rec_qual[irec] = qual_val;
299 }
300
301 }
302#else
303 double flux_val=0, confidence=0 ;
304 check( flux_val = cpl_image_get_interpolated( sci_img, fx, fy,
305 profile, radius, profile, radius, &confidence));
306 /* If error or undefined or pixel within kernel radius of edge */
307 if (confidence < 1) {
308 //rec_flux[irec] = 0;
309 rec_flux[irec] = flux_val;
310 rec_errs[irec] = 1;
311 rec_qual[irec] = QFLAG_MISSING_DATA;
312
313 /* If flux interpolation ok */
314 } else {
315
316 /* Now propagate errors */
317 float err_val = 0;
318 check( err_val = cpl_image_get_interpolated( var_img, fx, fy,
319 profile2, radius, profile2, radius, &confidence));
320
321 /* Now propagate flags */
322 int qual_val = 0;
323 for(j=yfirst;j<ylast-1;j++) {
324 j_nx=j*nx;
325 for(i=xfirst;i<xlast-1;i++) {
326 qual_val |= qual[j_nx+i];
327 }
328 }
329
330 /* Save the interpolated values */
331 rec_flux[irec] = flux_val;
332 rec_errs[irec] = sqrt(err_val);
333 rec_qual[irec] = qual_val;
334 }
335 }
336#endif
337
338 /* This is to perform flux conservation via the Jacobian */
339 rec_flux[irec] *= mult;
340 rec_errs[irec] *= mult;
341 }
342
343 cleanup: return;
344}
345/*****************************************************************************/
346
347/*****************************************************************************/
368/*****************************************************************************/
369
370static inline void
372 xsh_rec_list *rec_list,
373 int idx,
374 xsh_wavesol *wavesol,
375 xsh_xs_3 *model_config,
377 xsh_dispersol_list *disp_list,
378 float slit_min,
379 float slit_max,
380 double lambda_min,
381 int skip_low,
382 int skip_up,
383 xsh_rectify_param *rectify_par,
384 double slit_shift,
385 cpl_frame *slit_shiftab_frame)
386{
387
388 int ns, nl;
389 int nslit, nlambda ;
390 float * slit = NULL ;
391 double * lambda = NULL ;
392 //float * flux1 = NULL ;
393 //float * errs1 = NULL ;
394 //int * qual1 = NULL ;
395 int order, last_slit, last_lambda;
396 double *x_data = NULL;
397 double *y_data = NULL;
398
399 double *slit_shift_data = NULL;
400 const char *slit_shifttab_name = NULL;
401 cpl_table *shift_tab = NULL;
402 int shift_size;
403 double *slit_shifttab_data = NULL;
404 double *wave_shifttab_data = NULL;
405 cpl_image* var_img=NULL;
406
407 /* Check input parameters */
408 XSH_ASSURE_NOT_NULL( pre_sci);
409 XSH_ASSURE_NOT_NULL( rec_list);
410 check( slit = xsh_rec_list_get_slit( rec_list, idx));
411 XSH_ASSURE_NOT_NULL( slit);
412 check( lambda = xsh_rec_list_get_lambda( rec_list, idx));
413 XSH_ASSURE_NOT_NULL( lambda);
414 //check( flux1 = xsh_rec_list_get_data1( rec_list, idx));
415 //check( errs1 = xsh_rec_list_get_errs1( rec_list, idx));
416 //check( qual1 = xsh_rec_list_get_qual1( rec_list, idx));
417
418 var_img=cpl_image_duplicate(xsh_pre_get_errs( pre_sci));
419 cpl_image_multiply(var_img,var_img);
420
421 /* Create kernel profile (first time only ) */
422 if ( rec_profile == NULL && rectify_par->rectif_radius != 0 ) {
423 check( rec_profile = cpl_vector_new( CPL_KERNEL_DEF_SAMPLES));
424
425 assure(rectify_par->rectif_radius > 0,CPL_ERROR_ILLEGAL_INPUT,
426 "rectify-radius must be positive");
427 check( cpl_vector_fill_kernel_profile( rec_profile,
428 rectify_par->kernel_type, rectify_par->rectif_radius));
429 err_profile=cpl_vector_duplicate(rec_profile);
430 cpl_vector_multiply(err_profile,err_profile);
431 }
432
433 /* Get the order nb */
434 check( order = xsh_rec_list_get_order( rec_list, idx));
435 check( nslit = xsh_rec_list_get_nslit( rec_list, idx));
436 check( nlambda = xsh_rec_list_get_nlambda( rec_list, idx));
437
438 last_slit = nslit-1;
439 last_lambda = nlambda-1;
440
441 XSH_CALLOC( x_data, double, nslit*nlambda);
442 XSH_CALLOC( y_data, double, nslit*nlambda);
443
444 /* Fill slit and lambda arrays */
445 for( ns = 0 ; ns < nslit ; ns++ ) {
446 /* center the first pixel */
447 /*rec_list->list[idx].slit[ns] = slit_min +
448 ((float)ns*rectify_par->rectif_bin_space)+
449 rectify_par->rectif_bin_space/2.0;*/
450 rec_list->list[idx].slit[ns] = slit_min +
451 (double)ns*rectify_par->rectif_bin_space;
452 }
453
454 for( nl = 0 ; nl < nlambda ; nl++ ) {
455 /* center the first pixel */
456 rec_list->list[idx].lambda[nl] = lambda_min +
457 (double)nl*rectify_par->rectif_bin_lambda;
458 }
459
460 xsh_msg_dbg_low( "Order %d slit (%f,%f):%d lambda (%f,%f):%d",
461 order, rec_list->list[idx].slit[0],
462 rec_list->list[idx].slit[last_slit], nslit,
463 rec_list->list[idx].lambda[0],
464 rec_list->list[idx].lambda[last_lambda], nlambda);
465
466 /* Compute shift for wavelength */
467 XSH_CALLOC( slit_shift_data, double, nlambda);
468 if( slit_shiftab_frame != NULL){
469 check( slit_shifttab_name = cpl_frame_get_filename( slit_shiftab_frame));
470 xsh_msg_dbg_low("Load slit shift tab %s", slit_shifttab_name);
471 XSH_TABLE_LOAD( shift_tab, slit_shifttab_name);
472 check( wave_shifttab_data = cpl_table_get_data_double( shift_tab,
474 check( slit_shifttab_data = cpl_table_get_data_double( shift_tab,
476 shift_size = cpl_table_get_nrow( shift_tab);
477
478 for( nl = 0 ; nl < nlambda ; nl++ ) {
479 double wave = lambda[nl];
480 double cshift;
481 cshift = xsh_data_interpolate( wave, shift_size, wave_shifttab_data,
482 slit_shifttab_data);
483 slit_shift_data[nl] = cshift;
484 }
485
487 FILE *shift_file = NULL;
488 char shift_name[256];
489
490 sprintf( shift_name, "shift_%d_%.1f_%.1f.dat", order, slit_min, slit_max);
491 shift_file = fopen( shift_name, "w+");
492
493 fprintf( shift_file, "# wave delta_s\n");
494
495 for( nl = 0 ; nl < nlambda ; nl++ ) {
496 fprintf ( shift_file, "%f %f\n", lambda[nl], slit_shift_data[nl]);
497 }
498
499 fclose( shift_file);
500 }
501 }
502
503 /* GET X/Y positions */
504 if ( wavesol != NULL){
505 for( ns = skip_low ; ns < (nslit - skip_up) ; ns++ ) {
506 double xslit = slit[ns]+slit_shift;
507
508 for( nl = 0 ; nl < nlambda ; nl++ ) {
509 int irec;
510 double xlambda = lambda[nl];
511 double xslitrec;
512
513 irec= nl+ns*nlambda;
514 xslitrec = xslit+slit_shift_data[nl];
515 /* these functions make the conversion according to binning */
516 check( x_data[irec] = xsh_wavesol_eval_polx( wavesol, xlambda, order,
517 xslitrec));
518 check( y_data[irec] = xsh_wavesol_eval_poly( wavesol, xlambda, order,
519 xslitrec));
520 }
521 }
522 }
523 else{
524 XSH_ASSURE_NOT_NULL( model_config);
525 for( ns = skip_low ; ns < (nslit - skip_up) ; ns++ ) {
526 double xslit = slit[ns]+slit_shift;
527 int ns_nl=ns*nlambda;
528 for( nl = 0 ; nl < nlambda ; nl++ ) {
529 int irec;
530 double pm_x=0, pm_y=0;
531 double xlambda = lambda[nl];
532 double xslitrec;
533
534 irec= nl+ns_nl;
535 xslitrec = xslit+slit_shift_data[nl];
536 check( xsh_model_get_xy( model_config, instrument, xlambda, order,
537 xslitrec,
538 &pm_x, &pm_y));
539 x_data[irec] = convert_data_to_bin( pm_x, instrument->binx );
540 y_data[irec] = convert_data_to_bin( pm_y, instrument->biny );
541 }
542 }
543 }
544
546 FILE *rec_file = NULL;
547 char recfile_name[256];
548
549 sprintf( recfile_name, "rec_%d_%.1f_%.1f.dat", order, slit_min ,slit_max);
550 rec_file = fopen( recfile_name, "w+");
551
552 fprintf( rec_file, "#irec wave slit x y\n");
553
554 for( ns = skip_low ; ns < (nslit - skip_up) ; ns++ ) {
555 double xslit = slit[ns]+slit_shift;
556
557 for( nl = 0 ; nl < nlambda ; nl++ ) {
558 int irec;
559 double xlambda = lambda[nl];
560 double xslitrec;
561
562 irec= nl+ns*nlambda;
563 xslitrec = xslit+slit_shift_data[nl];
564 fprintf ( rec_file, "%d %f %f %f %f\n", irec, xlambda, xslitrec, x_data[irec],
565 y_data[irec]);
566 }
567 }
568 fclose( rec_file);
569 }
570
571 if (disp_list == NULL) {
572
573 for( ns = skip_low ; ns < (nslit - skip_up) ; ns++ ) {
574 int ns_nl=ns*nlambda;
575 for( nl = 0 ; nl < nlambda ; nl++ ) {
576 double fx=0, fy=0;
577 int irec=0;
578
579 irec= nl+ns_nl;
580 fx = x_data[irec];
581 fy = y_data[irec];
582
583 check( xsh_rec_list_rectify( rec_list, idx, irec, pre_sci,var_img,
584 fx, fy, rectify_par->rectif_radius,rec_profile,err_profile,1));
585 }
586
587 }
588 }
589 else{
590 cpl_polynomial *plambdadx = NULL, *plambdady = NULL;
591 cpl_polynomial *pslitdx = NULL, *pslitdy = NULL;
592 cpl_vector *val = NULL;
593 //int abs_order=0;
594
595 check( plambdadx = cpl_polynomial_duplicate(
596 disp_list->list[idx].lambda_poly));
597 check( plambdady = cpl_polynomial_duplicate(
598 disp_list->list[idx].lambda_poly));
599 check( pslitdx = cpl_polynomial_duplicate(
600 disp_list->list[idx].slit_poly));
601 check( pslitdy = cpl_polynomial_duplicate(
602 disp_list->list[idx].slit_poly));
603
604 //abs_order = disp_list->list[idx].absorder;
605
606 check( cpl_polynomial_derivative( plambdadx, 0));
607 check( cpl_polynomial_derivative( plambdady, 1));
608 check( cpl_polynomial_derivative( pslitdx, 0));
609 check( cpl_polynomial_derivative( pslitdy, 1));
610
611 check( val = cpl_vector_new(2));
612 double fx, fy;
613 double ldx, ldy, sdx, sdy, abs_determinant;
614 int irec;
615 int ns_nlambda=0;
616 double bin_size=rectify_par->rectif_bin_space*rectify_par->rectif_bin_lambda;
617 for( ns = skip_low ; ns < (nslit - skip_up) ; ns++ ) {
618 ns_nlambda=ns*nlambda;
619 for( nl = 0 ; nl < nlambda ; nl++ ) {
620
621 irec= nl+ns_nlambda;
622 fx = x_data[irec];
623 fy = y_data[irec];
624
625 check( cpl_vector_set( val, 0, fx));
626 check( cpl_vector_set( val, 1, fy));
627
628 check( ldx = cpl_polynomial_eval( plambdadx, val));
629 check( ldy = cpl_polynomial_eval( plambdady, val));
630 check( sdx = cpl_polynomial_eval( pslitdx, val));
631 check( sdy = cpl_polynomial_eval( pslitdy, val));
632 abs_determinant = bin_size/(fabs(ldx*sdy-sdx*ldy));
633
634 check( xsh_rec_list_rectify( rec_list, idx, irec, pre_sci,var_img,
635 fx, fy,
636 rectify_par->rectif_radius,
637 rec_profile, err_profile,abs_determinant));
638 }
639 }
640 xsh_free_polynomial( &plambdadx);
641 xsh_free_polynomial( &plambdady);
642 xsh_free_polynomial( &pslitdx);
643 xsh_free_polynomial( &pslitdy);
644 xsh_free_vector( &val);
645 }
646
647 cleanup:
648 xsh_free_image(&var_img);
651 XSH_FREE( x_data);
652 XSH_FREE( y_data);
653 XSH_FREE( slit_shift_data);
654 XSH_TABLE_FREE( shift_tab);
655 return ;
656}
657/*****************************************************************************/
658
659/*****************************************************************************/
660/*****************************************************************************/
670void
671xsh_get_slit_edges( cpl_frame *slitmap_frame, double *sdown, double *sup,
672 double *sldown, double *slup, xsh_instrument *instrument)
673{
674 const char *slitmap_name = NULL;
675 cpl_propertylist *slitmap_header = NULL;
676
677 if ( slitmap_frame != NULL) {
678 XSH_ASSURE_NOT_NULL( sdown);
680
681 check( slitmap_name = cpl_frame_get_filename( slitmap_frame));
682 check( slitmap_header = cpl_propertylist_load( slitmap_name, 0));
683
684 *sdown = xsh_pfits_get_slitmap_median_edglo( slitmap_header);
685 if ( cpl_error_get_code() != CPL_ERROR_NONE){
687 "Keyword 'MEDIAN EDGLO SLIT' not found in SLIT_MAP %s. Using default value %f",
688 slitmap_name, MIN_SLIT);
690 }
691
692 *sup = xsh_pfits_get_slitmap_median_edgup( slitmap_header);
693 if ( cpl_error_get_code() != CPL_ERROR_NONE){
695 "Keyword 'MEDIAN EDGUP SLIT' not found in SLIT_MAP %s. Using default value %f",
696 slitmap_name, MAX_SLIT);
698 }
700 XSH_ASSURE_NOT_NULL( sldown);
701 XSH_ASSURE_NOT_NULL( slup);
702 *sldown = xsh_pfits_get_slitmap_median_sliclo( slitmap_header);
703 if ( cpl_error_get_code() != CPL_ERROR_NONE){
705 "Keyword 'MEDIAN SLICLO SLIT' not found in SLIT_MAP %s. Using default value %f",
706 slitmap_name, SLITLET_CEN_CENTER-2);
708 }
709 *slup = xsh_pfits_get_slitmap_median_slicup( slitmap_header);
710 if ( cpl_error_get_code() != CPL_ERROR_NONE){
712 "Keyword 'MEDIAN SLICUP SLIT' not found in SLIT_MAP %s. Using default value %f",
713 slitmap_name, SLITLET_CEN_CENTER+2);
715 }
716 }
717 }
718 else {
719 xsh_msg_warning( "No provided SLIT_MAP. Using default values: 'MEDIAN EDGLO SLIT' %f - 'MEDIAN EDGUP SLIT' %f",
721 *sdown = MIN_SLIT;
722 *sup = MAX_SLIT;
724 xsh_msg_warning( "Using default values: 'MEDIAN SLICLO SLIT' %f - 'MEDIAN SLICUP SLIT' %f",
727 XSH_ASSURE_NOT_NULL( sldown);
728 XSH_ASSURE_NOT_NULL( slup);
729 *sldown = SLITLET_CEN_CENTER-2;
730 *slup = SLITLET_CEN_CENTER+2;
731 }
732 }
733
735 xsh_msg( "IFU limits: slitlet DOWN [%f %f], size: %f arcsec",
736 *sdown, *sldown, *sldown-*sdown);
737 xsh_msg( "IFU limits: slitlet CEN [%f %f], size: %f arcsec",
738 *sldown, *slup, *slup-*sldown);
739 xsh_msg( "IFU limits: slitlet UP [%f %f], size: %f arcsec",
740 *slup, *sup, *sup-*slup);
741 }
742 else{
743 xsh_msg( "SLIT limits: [%f %f], total size: %f arcsec", *sdown, *sup, *sup-*sdown);
744 }
745 cleanup:
746 xsh_free_propertylist( &slitmap_header);
747 return ;
748}
749/*****************************************************************************/
750
759 double *slit_min, int *nslit, XSH_MODE mode)
760{
761 double smin= MIN_SLIT, smax=MAX_SLIT;
762 double slit_step;
763
764 XSH_ASSURE_NOT_NULL( rectify_par);
765 XSH_ASSURE_NOT_NULL( slit_min);
766 XSH_ASSURE_NOT_NULL( nslit);
767
768 slit_step = rectify_par->rectif_bin_space;
769
770 if ( mode == XSH_MODE_SLIT){
771 if ( rectify_par->rectify_full_slit != 1 ) {
772 xsh_msg_warning(" Option not READY go to FULL_SLIT");
773 }
774 *nslit = (smax - smin)/slit_step;
775 *slit_min = smin;
776 smax = smin+(*nslit-1)*slit_step;
777
778 xsh_msg( "SLIT : (%.3f,%.3f) used only (%.3f,%.3f) in %d elts",
779 MIN_SLIT, MAX_SLIT, smin, smax, *nslit);
780 }
781
782 cleanup:
783 return;
784}
785/*****************************************************************************/
786
787static void
789 xsh_rectify_param * rectify_par )
790{
791 int i ;
792 double step = rectify_par->rectif_bin_lambda;
793
794 xsh_msg_dbg_high( "Adjust Lambdas: %d orders", spec_list->size ) ;
795
796 /* Adjust min/max according to lambda step */
797 for ( i = 0 ; i<spec_list->size ; i++ ) {
798 double min, max ;
799 int imin, imax ;
800
801 min = spec_list->list[i].lambda_min ;
802 max = spec_list->list[i].lambda_max ;
803 imin = ceil( min/step) ;
804 min = (double)imin * step;
805 imax = floor( max / step);
806 max = (double)imax*step ;
807 xsh_msg_dbg_high( "--- test lambdamin - %d : %.4lf -> %.4lf, %lf -> %lf", i,
808 spec_list->list[i].lambda_min, min,
809 spec_list->list[i].lambda_max, max ) ;
810 spec_list->list[i].lambda_min = min;
811 spec_list->list[i].lambda_max = max;
812 }
813}
814
815/* AMO: function not used
816static float
817get_step_slit( float * slit, int nslit )
818{
819 float result = 0. ;
820 float min, max ;
821
822 max = *(slit+nslit-1) ;
823 min = *slit ;
824 result = (max-min)/(float)nslit ;
825
826 xsh_msg( " Step Slit = %f", result ) ;
827 return result ;
828
829}
830*/
831
832/*****************************************************************************/
865/*****************************************************************************/
866cpl_frame *
867xsh_rectify( cpl_frame *sci_frame,
868 cpl_frame * orderlist_frame,
869 cpl_frame *wavesol_frame,
870 cpl_frame * model_frame,
872 xsh_rectify_param *rectify_par,
873 cpl_frame *spectralformat_frame,
874 cpl_frame *disp_tab_frame,
875 const char * res_name,
876 cpl_frame** res_frame_ext,
877 cpl_frame** res_frame_tab,
878 const char* rec_prefix)
879{
880 cpl_frame *result = NULL;
881 xsh_order_list *orderlist = NULL ;
882 int nslit;
883 double slit_min;
884 char tag[256];
885
886 XSH_ASSURE_NOT_NULL( orderlist_frame);
887
888 check( orderlist = xsh_order_list_load ( orderlist_frame, instrument));
889
890 sprintf(tag,"%s_%s",rec_prefix,
892
893 xsh_rec_slit_size( rectify_par, &slit_min, &nslit, XSH_MODE_SLIT);
894
895 check( result = xsh_rectify_orders( sci_frame, orderlist,
896 wavesol_frame, model_frame, instrument, rectify_par, spectralformat_frame,
897 disp_tab_frame, res_name, tag,
898 res_frame_ext, res_frame_tab,
899 0, 100, slit_min, nslit,0, NULL));
900
901 cleanup:
902 xsh_order_list_free( &orderlist);
903 return result;
904}
905/*****************************************************************************/
906
928cpl_frame * xsh_rectify_and_shift( cpl_frame *sci_frame,
929 cpl_frame * orderlist_frame,
930 cpl_frame *wavesol_frame,
931 cpl_frame * model_frame,
933 xsh_rectify_param *rectify_par,
934 cpl_frame *spectralformat_frame,
935 cpl_frame * loc_frame,
936 cpl_frame * loc0_frame,
937 double *throw_shift,
938 cpl_frame *disp_tab_frame,
939 const char * res_name,
940 cpl_frame** res_frame_ext,
941 cpl_frame** res_frame_tab)
942{
943 cpl_frame *result = NULL, *shift_rec = NULL;
944 cpl_frame * shift_rec_ext = NULL;
945 xsh_order_list *orderlist = NULL ;
946 int nslit;
947 double slit_min;
948 //const char *filename = NULL;
949 const char *tag = NULL;
950 cpl_propertylist *header = NULL;
951 double bin_space;
952 double shift_arcsec=0, shift_pix=0, shift_bin_arcsec=0;
953 double rectify_shift;
954#if 0
955 /* For test purpose */
956 cpl_frame * result0 = NULL;
957 char * res0_name = NULL ;
958#endif
959 xsh_msg( "===> xsh_rectify_and_shift" ) ;
960
961 XSH_ASSURE_NOT_NULL( orderlist_frame);
962 XSH_ASSURE_NOT_NULL( rectify_par);
963
964 check( orderlist = xsh_order_list_load ( orderlist_frame, instrument));
965
967 xsh_rec_slit_size( rectify_par, &slit_min, &nslit, XSH_MODE_SLIT);
968 bin_space = rectify_par->rectif_bin_space;
969
970 if ( throw_shift == NULL) {
971 if ( loc0_frame != NULL){
972 XSH_ASSURE_NOT_NULL( loc_frame ) ;
973 xsh_msg("Compute shift with localization");
974 shift_arcsec = compute_shift_with_localization( loc0_frame, loc_frame);
975 }
976 else{
977 xsh_msg("Fixed shift with localization to 0");
978 shift_arcsec = 0;
979 }
980 }
981 else {
982 shift_arcsec = *throw_shift;
983 xsh_msg("Use throw shift %f", shift_arcsec);
984 }
985
986 shift_pix = xsh_round_double(shift_arcsec/bin_space);
987 shift_bin_arcsec = shift_pix*bin_space;
988 rectify_shift = shift_bin_arcsec-shift_arcsec;
989
990 xsh_msg(" Mesured Shift for rectify : %f", rectify_shift);
991
992 check( shift_rec = xsh_rectify_orders( sci_frame, orderlist,
993 wavesol_frame, model_frame, instrument,
994 rectify_par, spectralformat_frame,
995 disp_tab_frame, res_name, tag,
996 &shift_rec_ext, res_frame_tab,
997 0, 100, slit_min, nslit, rectify_shift,
998 NULL));
999
1000 check( result = xsh_shift( shift_rec, instrument, res_name, shift_bin_arcsec,
1001 res_frame_ext));
1002
1003 cleanup:
1004 xsh_free_frame( &shift_rec);
1005 xsh_free_frame( &shift_rec_ext);
1006 xsh_order_list_free( &orderlist);
1007 xsh_free_propertylist( &header);
1008 return result;
1009}
1010/*****************************************************************************/
1011
1012/*****************************************************************************/
1013/*****************************************************************************/
1014cpl_frameset* xsh_rectify_ifu(cpl_frame *sci_frame,
1015 cpl_frame *orderlist_frame, cpl_frameset *wavesol_frameset,
1016 cpl_frameset *shiftifu_frameset,
1017 cpl_frame *model_frame, xsh_instrument *instrument,
1018 xsh_rectify_param *rectify_par, cpl_frame *spectralformat_frame,
1019 cpl_frame * slitmap_frame,
1020 cpl_frameset** rec_frameset_ext,
1021 cpl_frameset** rec_frameset_tab,
1022 const char* rec_prefix )
1023{
1024 cpl_frameset *result = NULL;
1025 xsh_order_list *orderlist = NULL ;
1026
1027 XSH_ASSURE_NOT_NULL( orderlist_frame);
1028 check( orderlist = xsh_order_list_load ( orderlist_frame, instrument));
1029
1030
1031 XSH_REGDEBUG("TODO : ADD disp_tab frameset, res_frame_ext frameset");
1032
1033 check( result = xsh_rectify_orders_ifu( sci_frame, orderlist,
1034 wavesol_frameset, shiftifu_frameset, model_frame, instrument, rectify_par,
1035 spectralformat_frame,
1036 slitmap_frame,rec_frameset_ext,rec_frameset_tab, 0, 100, rec_prefix));
1037
1038 cleanup:
1039 xsh_order_list_free( &orderlist);
1040 return result;
1041}
1042/*****************************************************************************/
1043
1044/*****************************************************************************/
1076/*****************************************************************************/
1077/* TODO: clarify this function */
1078cpl_frame*
1079xsh_rectify_orders( cpl_frame *sci_frame,
1080 xsh_order_list *orderlist,
1081 cpl_frame *wavesol_frame,
1082 cpl_frame *model_frame,
1084 xsh_rectify_param *rectify_par,
1085 cpl_frame *spectralformat_frame,
1086 cpl_frame *disp_tab_frame,
1087 const char *res_name,
1088 const char* tag,
1089 cpl_frame** res_frame_ext,
1090 cpl_frame** res_frame_tab,
1091 int min_index,
1092 int max_index,
1093 double slit_min,
1094 int nslit,
1095 double slit_shift,
1096 cpl_frame *slitshift_tab)
1097{
1098 xsh_rec_list * rec_list = NULL ;
1099 xsh_xs_3 config_model ;
1100 //xsh_xs_3* p_xs_3;
1101 int solution_type = 0;
1102 xsh_wavesol * wavesol = NULL ;
1103 cpl_frame * res_frame = NULL ;
1104 xsh_pre * pre_sci = NULL ;
1105 int i, nlambda, order ;
1106 float slit_max ;
1107 xsh_spectralformat_list *spec_list = NULL ;
1108 char tag_drl[256];
1109 char res_name_drl[256];
1110 //cpl_mask* qual_mask = NULL;
1111 xsh_dispersol_list* disp_list = NULL;
1112 int rec_min_index, rec_max_index;
1113 const char* solution_type_name[2] = { "POLY", "MODEL"};
1114 int found_line=true;
1115 double barycor=0;
1116 double helicor=0;
1117
1118 XSH_ASSURE_NOT_NULL( sci_frame);
1119 XSH_ASSURE_NOT_NULL( orderlist);
1120 XSH_ASSURE_NOT_NULL( spectralformat_frame);
1121 XSH_ASSURE_NOT_NULL( rectify_par);
1123 XSH_ASSURE_NOT_ILLEGAL( min_index <= max_index);
1124
1125 if (rectify_par->conserve_flux){
1126 assure( disp_tab_frame != NULL,CPL_ERROR_ILLEGAL_INPUT,
1127 "If -rectify-conserve-flux=TRUE "
1128 "you must provide an input dispersion table "
1129 "tagged as DISP_TAB_%s or DISP_TAB_AFC_%s .",
1132 check( disp_list = xsh_dispersol_list_load( disp_tab_frame, instrument));
1133 }
1134
1135 check( pre_sci = xsh_pre_load( sci_frame, instrument));
1136 check( rec_list = xsh_rec_list_create( instrument));
1137
1138 check( xsh_baryvel(pre_sci->data_header, &barycor,&helicor));
1139 cpl_propertylist_append_double(pre_sci->data_header,XSH_QC_VRAD_BARYCOR,barycor);
1140 cpl_propertylist_set_comment(pre_sci->data_header,XSH_QC_VRAD_BARYCOR,
1142 cpl_propertylist_append_double(pre_sci->data_header,XSH_QC_VRAD_HELICOR,helicor);
1143 cpl_propertylist_set_comment(pre_sci->data_header,XSH_QC_VRAD_HELICOR,
1145
1146 if ( model_frame != NULL) {
1147 solution_type = XSH_RECTIFY_TYPE_MODEL;
1148 check(xsh_model_temperature_update_frame(&model_frame,sci_frame,
1149 instrument,&found_line));
1150 check( xsh_model_config_load_best( model_frame, &config_model));
1151 }
1152 else if ( wavesol_frame != NULL ) {
1153 solution_type = XSH_RECTIFY_TYPE_POLY;
1154 check( wavesol = xsh_wavesol_load( wavesol_frame, instrument));
1155 }
1156 else {
1158 "Undefined solution type (POLY or MODEL). See your input sof");
1159 }
1160
1161 xsh_msg("===== Solution type %s", solution_type_name[solution_type]);
1162
1163 check( spec_list = xsh_spectralformat_list_load( spectralformat_frame,
1164 instrument));
1165
1166 /* Parameters */
1167 xsh_msg_dbg_high( "Kernel: %s, Radius: %lf", rectify_par->rectif_kernel,
1168 rectify_par->rectif_radius);
1169 xsh_msg_dbg_high( "Slit Binning: %f, Lambda Binning: %f",
1170 rectify_par->rectif_bin_space,
1171 rectify_par->rectif_bin_lambda);
1172
1173 /* Populate list */
1174 adjust_lambdas( spec_list, rectify_par);
1175
1176 /* Should make sure that we have the same nb of slits */
1177 rec_list->nslit = nslit;
1178 rec_list->slit_min = slit_min;
1179 slit_max = slit_min+rectify_par->rectif_bin_space*nslit;
1180 rec_list->slit_max = slit_max;
1181
1182 if (min_index >= 0){
1183 rec_min_index = min_index;
1184 }
1185 else{
1186 rec_min_index = 0;
1187 }
1188 if ( max_index < rec_list->size){
1189 rec_max_index = max_index;
1190 }
1191 else{
1192 rec_max_index = rec_list->size-1;
1193 }
1194
1195 /* AMO added to make sure bad pixels are flagged before re-sampling */
1196
1197 //xsh_msg("ok1");
1198 //int decode_bp=2147483647;
1199 int decode_bp=instrument->decode_bp;
1200 cpl_mask* mask = xsh_qual_to_cpl_mask(pre_sci->qual,decode_bp);
1201 cpl_image_reject_from_mask(pre_sci->data,mask);
1202 //cpl_detector_interpolate_rejected(pre_sci->data);
1203 xsh_free_mask(&mask);
1204
1205
1206 for( i = rec_min_index; i <= rec_max_index; i++ ) {
1207 double lambda_min, lambda_max ;
1208 double n;
1209
1210 order = spec_list->list[i].absorder ;
1211 rec_list->list[i].nslit = nslit ;
1212 lambda_min = spec_list->list[i].lambda_min ;
1213 lambda_max = spec_list->list[i].lambda_max ;
1214
1215 /* Calculate nlambda using lambda min/max and rectif_lambda-space */
1216 n = (lambda_max-lambda_min)/rectify_par->rectif_bin_lambda+1.0;
1217 nlambda = (int)xsh_round_double(n);
1218 assure(nlambda>0,CPL_ERROR_ILLEGAL_INPUT,
1219 "Number of wavelength sampling points is %d. Must be > 0. "
1220 "You may have to decrease rectify-bin-lambda param value.",
1221 nlambda);
1222
1223 assure(nslit>0,CPL_ERROR_ILLEGAL_INPUT,
1224 "Number of spatial sampling points is %d. Must be > 0. "
1225 "You may have to decrease rectify-bin-slit param value.",
1226 nslit);
1227
1228 rec_list->list[i].nlambda = nlambda ;
1229
1230
1231 check( xsh_rec_list_set_data_size( rec_list, i, order, nlambda, nslit));
1232 /*
1233 For each lambda/slit calculate X and Y and fill the grid
1234 with flux, errs, qual
1235 Needs wavesol
1236 */
1237 xsh_msg_dbg_medium( "Order %d, Nslit: %d, Nlambda: %d", order, nslit, nlambda ) ;
1238 xsh_msg_dbg_medium( " Lambda (%f,%f) Slit(%f,%f)",
1239 lambda_min, lambda_max , slit_min, slit_max);
1240
1241 //check( qual_mask = xsh_pre_get_bpmap( pre_sci));
1242 /* in case of simple data re-sampling we do not need to declare bad pixels in the image data extension
1243 * because we simply re-map from one space to another
1244 check( cpl_image_reject_from_mask( pre_sci->data, qual_mask));
1245 */
1246 check( fill_rectified( pre_sci, rec_list, i, wavesol,
1247 &config_model, instrument, disp_list,
1248 slit_min, slit_max, lambda_min, 0, 0, rectify_par,
1249 slit_shift, slitshift_tab));
1250 }
1251
1252 rec_list->slit_min = rec_list->list[rec_min_index].slit[0];
1253 rec_list->slit_max = rec_list->list[rec_min_index].slit[nslit-1];
1254
1255 xsh_msg_dbg_medium( "Saving Rectified Frame '%s'", res_name);
1256
1257 sprintf( tag_drl ,"%s_DRL", tag);
1258 sprintf( res_name_drl ,"DRL_%s", res_name);
1259
1260
1261 check( xsh_rec_list_update_header( rec_list, pre_sci, rectify_par, tag_drl));
1262
1263 if ( slitshift_tab != NULL){
1264 xsh_frame_set_shiftifu_ref( rec_list->header, slitshift_tab);
1265 }
1266
1267 check( res_frame = xsh_rec_list_save( rec_list, res_name_drl, tag_drl, CPL_TRUE));
1268 xsh_add_temporary_file(res_name_drl);
1269
1270 if ( res_frame_ext != NULL){
1271
1272 check( *res_frame_ext = xsh_rec_list_save2(rec_list,res_name,tag));
1273 sprintf( tag_drl ,"%s_TAB", tag);
1274 sprintf( res_name_drl ,"TAB_%s", res_name);
1275 check( *res_frame_tab = xsh_rec_list_save_table(rec_list,res_name_drl,tag_drl,
1276 CPL_TRUE));
1277
1278 xsh_add_temporary_file(res_name_drl);
1279 }
1280
1281 cleanup :
1282 xsh_dispersol_list_free( &disp_list);
1283 xsh_rec_list_free( &rec_list);
1284 xsh_wavesol_free( &wavesol);
1285 xsh_spectralformat_list_free( &spec_list);
1286 xsh_pre_free( &pre_sci);
1287 return res_frame;
1288}
1289/*****************************************************************************/
1290
1291static void xsh_get_shift_ref( cpl_frameset *set,
1292 double *down, double *cen, double *up)
1293{
1294 cpl_frame *offsetdown_frame = NULL;
1295 cpl_frame *offsetcen_frame = NULL;
1296 cpl_frame *offsetup_frame = NULL;
1297 const char *offsetdown_name = NULL;
1298 const char *offsetcen_name = NULL;
1299 const char *offsetup_name = NULL;
1300 cpl_propertylist *offsetdown_list = NULL;
1301 cpl_propertylist *offsetcen_list = NULL;
1302 cpl_propertylist *offsetup_list = NULL;
1303
1304 XSH_ASSURE_NOT_NULL( down);
1305 XSH_ASSURE_NOT_NULL( cen);
1307
1308 check( offsetdown_frame = cpl_frameset_get_frame( set, 0));
1309 check( offsetcen_frame = cpl_frameset_get_frame( set, 1));
1310 check( offsetup_frame = cpl_frameset_get_frame( set, 2));
1311 check( offsetdown_name = cpl_frame_get_filename( offsetdown_frame));
1312 check( offsetcen_name = cpl_frame_get_filename( offsetcen_frame));
1313 check( offsetup_name = cpl_frame_get_filename( offsetup_frame));
1314 check( offsetdown_list = cpl_propertylist_load( offsetdown_name, 0));
1315 check( offsetcen_list = cpl_propertylist_load( offsetcen_name, 0));
1316 check( offsetup_list = cpl_propertylist_load( offsetup_name, 0));
1317
1318 check( *down = xsh_pfits_get_shiftifu_slitref( offsetdown_list));
1319 check( *cen = xsh_pfits_get_shiftifu_slitref( offsetcen_list));
1320 check( *up = xsh_pfits_get_shiftifu_slitref( offsetup_list));
1321
1322 cleanup:
1323 xsh_free_propertylist( &offsetdown_list);
1324 xsh_free_propertylist( &offsetcen_list);
1325 xsh_free_propertylist( &offsetup_list);
1326 return;
1327}
1328/*****************************************************************************/
1329
1330
1331/*****************************************************************************/
1332static void xsh_frame_set_shiftifu_ref( cpl_propertylist *header,
1333 cpl_frame *shift_frame)
1334{
1335 cpl_propertylist *shift_header = NULL;
1336 const char *shift_name = NULL;
1337 double waveref, slitref;
1338
1339 XSH_ASSURE_NOT_NULL( header);
1340 XSH_ASSURE_NOT_NULL( shift_frame);
1341
1342 check( shift_name = cpl_frame_get_filename( shift_frame));
1343 check( shift_header = cpl_propertylist_load( shift_name, 0));
1344
1345 check( waveref = xsh_pfits_get_shiftifu_lambdaref( shift_header));
1346 check( slitref = xsh_pfits_get_shiftifu_slitref( shift_header));
1347
1348 check( xsh_pfits_set_shiftifu_lambdaref( header, waveref));
1349 check( xsh_pfits_set_shiftifu_slitref( header, slitref));
1350
1351 cleanup:
1352 xsh_free_propertylist( &shift_header);
1353 return;
1354}
1355/*****************************************************************************/
1356
1357/*****************************************************************************/
1372void xsh_compute_slitlet_limits( cpl_frameset *shift_set, double sdown,
1373 double sldown, double slup, double sup, double slit_bin,
1374 double *slitmin_tab, int *nslit_tab, double *slitcen_tab)
1375{
1376 double offset_down=0, offset_cen=0, offset_up=0;
1377 double cen_slitlet_down_size;
1378 double cen_slitlet_up_size=0;
1379 double down_slitlet_down_size;
1380 double down_slitlet_up_size;
1381 double up_slitlet_down_size;
1382 double up_slitlet_up_size=0;
1383 double slitlet_up_size;
1384 double slitlet_down_size;
1385 double cen_ref_down, cen_ref_up;
1386 double cen_pixpos = 0.0;
1387 int ref_slit_size;
1388 int down_nslit;
1389 int up_nslit;
1390
1391 XSH_ASSURE_NOT_NULL( shift_set);
1392 XSH_ASSURE_NOT_NULL( slitmin_tab);
1393 XSH_ASSURE_NOT_NULL( nslit_tab);
1394 XSH_ASSURE_NOT_NULL( slitcen_tab);
1395
1396 check( xsh_get_shift_ref( shift_set, &offset_down, &offset_cen,
1397 &offset_up));
1398
1399
1400 xsh_msg( "Shift reference: down %f arcsec, center %f arcsec, up %f arcsec",
1401 offset_down, offset_cen, offset_up);
1402
1403 slitcen_tab[0] = offset_down;
1404 slitcen_tab[1] = offset_cen;
1405 slitcen_tab[2] =offset_up;
1406
1407 down_slitlet_down_size = offset_down-sdown;
1408 down_slitlet_up_size = sldown-offset_down;
1409
1410 xsh_msg_dbg_medium("In down slitlet [%f,%f] size lo %f up %f", sdown, sldown,
1411 down_slitlet_down_size, down_slitlet_up_size);
1412
1413 cen_slitlet_down_size = offset_cen-sldown;
1414 cen_slitlet_up_size = slup-offset_cen;
1415
1416 xsh_msg_dbg_medium("In cen slitlet [%f,%f] size lo %f up %f", sldown, slup,
1417 cen_slitlet_down_size, cen_slitlet_up_size);
1418
1419 up_slitlet_down_size = offset_up-slup;
1420 up_slitlet_up_size = sup-offset_up;
1421
1422 xsh_msg_dbg_medium("In up slitlet [%f,%f] size lo %f up %f", slup, sup,
1423 up_slitlet_down_size, up_slitlet_up_size);
1424
1425 /* Slitlet down size : minimum between DOWN of center sliltet and UP of DOWN and up slitlets */
1426 if ( cen_slitlet_down_size < up_slitlet_up_size){
1427 slitlet_down_size = cen_slitlet_down_size;
1428 }
1429 else{
1430 slitlet_down_size = up_slitlet_up_size;
1431 }
1432 if ( slitlet_down_size > down_slitlet_up_size){
1433 slitlet_down_size = down_slitlet_up_size;
1434 }
1435 /* Slitlet up size : minimum between UP of center sliltet and DOWN of DOWN and up slitlets */
1436 if ( cen_slitlet_up_size < up_slitlet_down_size){
1437 slitlet_up_size = cen_slitlet_up_size;
1438 }
1439 else{
1440 slitlet_up_size = up_slitlet_down_size;
1441 }
1442 if ( slitlet_up_size > down_slitlet_down_size){
1443 slitlet_up_size = down_slitlet_down_size;
1444 }
1445
1446 xsh_msg_dbg_medium( "Slitlet size DOWN %f UP %f", slitlet_down_size, slitlet_up_size);
1447
1448 /* New version */
1449 cen_ref_down = offset_cen-slitlet_down_size;
1450 cen_ref_up = offset_cen+slitlet_up_size;
1451 if (cen_ref_down > 0){
1452 down_nslit = ceil( cen_ref_down/slit_bin);
1453 }
1454 else{
1455 down_nslit = floor( cen_ref_down/slit_bin);
1456 }
1457 if ( cen_ref_up > 0){
1458 up_nslit = ceil( cen_ref_up/slit_bin);
1459 }
1460 else{
1461 up_nslit = floor( cen_ref_up/slit_bin);
1462 }
1463
1464 xsh_msg( "Adjust central reference slitlet [%f %f] with bin %f arcsec: [%f %f]",
1465 cen_ref_down, cen_ref_up, slit_bin, down_nslit*slit_bin, up_nslit*slit_bin);
1466
1467 cen_pixpos = (offset_cen-down_nslit*slit_bin)/slit_bin;
1468 ref_slit_size = up_nslit-down_nslit+1;
1469 xsh_msg( "Center position in pixel %f (Total %d pix)", cen_pixpos, ref_slit_size);
1470 slitmin_tab[1] = down_nslit*slit_bin;
1471 nslit_tab[1] = ref_slit_size;
1472
1473 slitmin_tab[0] = offset_down-(ref_slit_size-1-cen_pixpos)*slit_bin;
1474 nslit_tab[0] = ref_slit_size;
1475
1476 slitmin_tab[2] = offset_up-(ref_slit_size-1-cen_pixpos)*slit_bin;
1477 nslit_tab[2] = ref_slit_size;
1478
1479 xsh_msg("Prepare the cube bin %f arcsec", slit_bin);
1480 xsh_msg("DOWN [%f, %f]", slitmin_tab[0], slitmin_tab[0]+nslit_tab[0]*slit_bin);
1481 xsh_msg("CEN [%f, %f]", slitmin_tab[1], slitmin_tab[1]+nslit_tab[1]*slit_bin);
1482 xsh_msg("UP [%f, %f]", slitmin_tab[2], slitmin_tab[2]+nslit_tab[2]*slit_bin);
1483
1484 cleanup:
1485 return;
1486}
1487/*****************************************************************************/
1488
1489/*****************************************************************************/
1512/*****************************************************************************/
1513cpl_frameset *
1514xsh_rectify_orders_ifu(cpl_frame *sci_frame,
1515 xsh_order_list *orderlist,
1516 cpl_frameset *wavesol_frameset,
1517 cpl_frameset *shift_frameset,
1518 cpl_frame *model_frame,
1520 xsh_rectify_param *rectify_par,
1521 cpl_frame *spectralformat_frame,
1522 cpl_frame * slitmap_frame,
1523 cpl_frameset ** res_frameset_ext,
1524 cpl_frameset ** res_frameset_tab,
1525 int min_index,
1526 int max_index,
1527 const char* rec_prefix)
1528{
1529 cpl_frameset *res_frameset = NULL ;
1530
1531 cpl_frame *wavesolup_frame = NULL;
1532 cpl_frame *wavesolcen_frame = NULL;
1533 cpl_frame *wavesoldown_frame = NULL;
1534 cpl_frame *shiftup_frame = NULL;
1535 cpl_frame *shiftcen_frame = NULL;
1536 cpl_frame *shiftdown_frame = NULL;
1537
1538 int slitlet=0;
1539 double sdown=0, sldown=0, slup=0, sup=0;
1540 //double slitlet_down_l=0;
1541 double slitlet_down_u=0;
1542
1543 double slitlet_cen_l=0, slitlet_cen_u=0;
1544 double slitlet_cen_u_bin=0;
1545 //double slitlet_up_u=0;
1546 double down_offset=0.0;
1547 int nslit=0;
1548 double bin_space;
1549 //XSH_ARM arm ;
1550 int nslit_tab[3];
1551 double slitmin_tab[3];
1552 double slitcen_tab[3];
1553 int cut_edges = 0;
1554
1555 XSH_ASSURE_NOT_NULL( sci_frame);
1556 XSH_ASSURE_NOT_NULL( orderlist);
1557 XSH_ASSURE_NOT_NULL( rectify_par);
1559 XSH_ASSURE_NOT_NULL( spectralformat_frame);
1560
1561 bin_space = rectify_par->rectif_bin_space;
1562 /*cut_edges = rectify_par->cut_edges;*/
1563
1564 //arm = xsh_instrument_get_arm( instrument ) ;
1565
1566 if ( wavesol_frameset != NULL){
1567 check( wavesoldown_frame = cpl_frameset_get_frame( wavesol_frameset, 0));
1568 check( wavesolcen_frame = cpl_frameset_get_frame( wavesol_frameset, 1));
1569 check( wavesolup_frame = cpl_frameset_get_frame( wavesol_frameset, 2));
1570 }
1571
1572 if ( shift_frameset != NULL){
1573 check( shiftdown_frame = cpl_frameset_get_frame( shift_frameset, 0));
1574 check( shiftcen_frame = cpl_frameset_get_frame( shift_frameset, 1));
1575 check( shiftup_frame = cpl_frameset_get_frame( shift_frameset, 2));
1576 }
1577
1578 /* Create result frameset */
1579 check( res_frameset = cpl_frameset_new());
1580
1581 check( xsh_get_slit_edges( slitmap_frame, &sdown,
1582 &sup, &sldown, &slup, instrument));
1583
1584 if ( shift_frameset != NULL){
1585 check( xsh_compute_slitlet_limits( shift_frameset, sdown,
1586 sldown, slup, sup, bin_space,
1587 slitmin_tab, nslit_tab, slitcen_tab));
1588 }
1589 else{
1590 slitlet_cen_l = sldown;
1591 slitlet_cen_u = slup;
1592
1593 //slitlet_down_l = sdown;
1594 slitlet_down_u = slitlet_cen_l;
1595 //slitlet_up_u = sup;
1596
1597 /* center */
1598 slitmin_tab[1] = slitlet_cen_l;
1599 nslit_tab[1] = ceil((slitlet_cen_u-slitlet_cen_l)/rectify_par->rectif_bin_space);
1600 slitlet_cen_u_bin = slitmin_tab[1]+nslit_tab[1]*rectify_par->rectif_bin_space;
1601
1602 /* down */
1603 nslit_tab[0] = nslit_tab[1];
1604 slitmin_tab[0] = slitlet_down_u-down_offset-
1605 nslit_tab[0]*rectify_par->rectif_bin_space;
1606
1607 /* up */
1608 slitmin_tab[2] = slitlet_cen_u_bin;
1609 nslit_tab[2] = nslit_tab[1];
1610 }
1611
1612 for ( slitlet = LOWER_IFU_SLITLET ; slitlet <= UPPER_IFU_SLITLET ;
1613 slitlet++ ) {
1614
1615 cpl_frame *wavesol_frame = NULL;
1616 cpl_frame *shift_frame = NULL;
1617 cpl_frame *res_frame = NULL;
1618 cpl_frame *resext_frame = NULL;
1619 cpl_frame *restab_frame = NULL;
1620 double slitlet_slit_min = 0.0;
1621 char tag[256];
1622 char res_name[256];
1623 double rec_shift=0.0;
1624 double rec_min=0.0;
1625 const char *slitlet_name = NULL;
1626
1627 switch( slitlet ) {
1628 case LOWER_IFU_SLITLET:
1629
1630 slitlet_name = "DOWN";
1631 nslit = nslit_tab[0];
1632 slitlet_slit_min = slitmin_tab[0];
1633 wavesol_frame = wavesoldown_frame;
1634 shift_frame = shiftdown_frame;
1635 sprintf(tag,"%s_%s",rec_prefix,
1636 XSH_GET_TAG_FROM_ARM(XSH_ORDER2D_DOWN_IFU,instrument));
1637 sprintf(res_name,"%s.fits",tag);
1638 break ;
1639
1640 case CENTER_IFU_SLITLET:
1641 slitlet_name = "CEN";
1642 slitlet_slit_min = slitmin_tab[1];
1643 nslit = nslit_tab[1];
1644 wavesol_frame = wavesolcen_frame;
1645 shift_frame = shiftcen_frame;
1646 sprintf(tag,"%s_%s",rec_prefix,
1647 XSH_GET_TAG_FROM_ARM(XSH_ORDER2D_CEN_IFU,instrument));
1648 sprintf(res_name,"%s.fits",tag);
1649 break ;
1650
1651 case UPPER_IFU_SLITLET:
1652 slitlet_name = "UP";
1653 slitlet_slit_min = slitmin_tab[2];
1654 nslit = nslit_tab[2];
1655
1656 wavesol_frame = wavesolup_frame;
1657 shift_frame = shiftup_frame;
1658 sprintf(tag,"%s_%s",rec_prefix,
1659 XSH_GET_TAG_FROM_ARM(XSH_ORDER2D_UP_IFU,instrument));
1660 sprintf(res_name,"%s.fits",tag);
1661 break ;
1662 }
1663 /* compute rec_shift */
1664 rec_min = slitlet_slit_min;
1665 rec_shift = 0;
1666
1667 xsh_msg("Cut edges with %d pix", cut_edges);
1668 rec_min += cut_edges*rectify_par->rectif_bin_space;
1669 nslit -= cut_edges*2;
1670
1671 XSH_CMP_INT( nslit, >, 0, "Check size in slit",);
1672
1673
1674 xsh_msg( "%s [%f %f] in %d pix of %f arcsec" ,slitlet_name,
1675 slitlet_slit_min,
1676 slitlet_slit_min+nslit*rectify_par->rectif_bin_space,
1677 nslit, rectify_par->rectif_bin_space);
1678
1679 check( res_frame = xsh_rectify_orders( sci_frame, orderlist,
1680 wavesol_frame, model_frame, instrument, rectify_par,
1681 spectralformat_frame, NULL, res_name, tag, &resext_frame, &restab_frame,
1682 min_index, max_index, rec_min, nslit, rec_shift,shift_frame));
1683
1684 check( cpl_frameset_insert( res_frameset, res_frame));
1685
1686 if(res_frameset_ext !=NULL && res_frameset_tab !=NULL) {
1687 check( cpl_frameset_insert(*res_frameset_ext, resext_frame));
1688 check( cpl_frameset_insert(*res_frameset_tab, restab_frame));
1689 }
1690 }
1691
1692 cleanup:
1693 return res_frameset ;
1694}
1695/*****************************************************************************/
1696
1697
1698/*****************************************************************************/
1706/*****************************************************************************/
1707static double
1709 cpl_frame *loc0_frame)
1710{
1711 double shift_a=0, shift_b=0, shift=0;
1712 xsh_localization *loc = NULL, *loc0 = NULL;
1713
1714 XSH_ASSURE_NOT_NULL( loc_frame);
1715 XSH_ASSURE_NOT_NULL( loc0_frame);
1716
1717 check( loc = xsh_localization_load( loc_frame));
1718 check( loc0 = xsh_localization_load( loc0_frame));
1719
1720 /* We assume here that the degree of the localization polynomial is 0 */
1721 check( shift_a = cpl_polynomial_eval_1d( loc0->cenpoly, 0., NULL));
1722 check( shift_b = cpl_polynomial_eval_1d( loc->cenpoly, 0., NULL));
1723
1724 xsh_msg_dbg_medium("Shift A %f B %f", shift_a, shift_b)
1725 shift = shift_b-shift_a;
1726
1727 xsh_msg( "Shift (localization) = %lf", shift);
1728
1729 cleanup:
1730 xsh_localization_free( &loc);
1731 xsh_localization_free( &loc0);
1732 return shift;
1733}
1734/*****************************************************************************/
1735/* AMO: function not used. For the moment commented out
1736static cpl_frame *
1737shift_with_localization( cpl_frame * rec_frame,
1738 cpl_frame * loc_frame,
1739 cpl_frame * loc0_frame,
1740 xsh_instrument * instrument,
1741 const char * fname,
1742 cpl_frame** res_frame_ext )
1743{
1744 cpl_frame * result = NULL ; //< The result frame (the
1745 //shifted B-A rectified image)
1746 xsh_rec_list * shifted_rec_list = NULL ; ///< The result list. local
1747 xsh_rec_list * rec_list = NULL;
1748 xsh_localization *loc_list = NULL;
1749 xsh_localization *loc0_list = NULL;
1750 float step_slit ; //< A-B object position ( ... )
1751 int nb_orders ; //< Nb of orders in B-A
1752 int order ;
1753 int absorder, nlambda, nslit ;
1754 float * flux1 = NULL ;
1755 float * errs1 = NULL ;
1756 int * qual1 = NULL ;
1757 float * shifted_flux1 = NULL ;
1758 float * shifted_errs1 = NULL ;
1759 int * shifted_qual1 = NULL ;
1760 float * slit = NULL ;
1761 double * lambda = NULL ;
1762 float * shifted_slit = NULL ;
1763 double * shifted_lambda = NULL ;
1764 const char* tag=NULL;
1765 const char* tag_drl=NULL;
1766 const char* fname_drl=NULL;
1767
1768 XSH_ASSURE_NOT_NULL( rec_frame ) ;
1769 XSH_ASSURE_NOT_NULL( loc_frame ) ;
1770 XSH_ASSURE_NOT_NULL( loc0_frame ) ;
1771 XSH_ASSURE_NOT_NULL( instrument ) ;
1772
1773 xsh_msg_dbg_low( "Entering shift_with_localization" ) ;
1774
1775 // Load needed frames into lists
1776 check( rec_list = xsh_rec_list_load( rec_frame, instrument));
1777 check( loc_list = xsh_localization_load( loc_frame));
1778 check( loc0_list = xsh_localization_load( loc0_frame));
1779 // Create shifted rectified list
1780 check( shifted_rec_list = xsh_rec_list_duplicate( rec_list, instrument ) ) ;
1781 xsh_msg( "Nb of orders in new rec list: %d", shifted_rec_list->size ) ;
1782 // Get the nb of orders
1783 check( nb_orders = rec_list->size ) ;
1784
1785 // Loop over orders
1786 for( order = 0 ; order < nb_orders ; order++ ) {
1787 int ns, nl ;
1788 int nshift, max_size ;
1789
1790 check( absorder = xsh_rec_list_get_order( rec_list, order ) ) ;
1791 check( nslit = xsh_rec_list_get_nslit( rec_list, order ) ) ;
1792 check( nlambda = xsh_rec_list_get_nlambda( rec_list, order ) ) ;
1793 xsh_msg( " Absolute order: %d, Nslit: %d, Nlambda: %d", absorder,
1794 nslit, nlambda ) ;
1795 max_size = nslit*nlambda ;
1796
1797 // Prepare shifting
1798 check( flux1 = xsh_rec_list_get_data1( rec_list, order ) ) ;
1799 check( qual1 = xsh_rec_list_get_qual1( rec_list, order ) ) ;
1800 check( errs1 = xsh_rec_list_get_errs1( rec_list, order ) ) ;
1801 check( slit = xsh_rec_list_get_slit( rec_list, order ) ) ;
1802 check( lambda = xsh_rec_list_get_lambda( rec_list, order ) ) ;
1803
1804 check( shifted_slit = xsh_rec_list_get_slit( shifted_rec_list,
1805 order ) ) ;
1806 check( shifted_lambda = xsh_rec_list_get_lambda( shifted_rec_list,
1807 order ) ) ;
1808 check( shifted_flux1 = xsh_rec_list_get_data1( shifted_rec_list,
1809 order ) ) ;
1810 check( shifted_qual1 = xsh_rec_list_get_qual1( shifted_rec_list,
1811 order ) ) ;
1812 check( shifted_errs1 = xsh_rec_list_get_errs1( shifted_rec_list,
1813 order ) ) ;
1814 // Set shifted list slits and lambda
1815 for( ns = 0 ; ns<nslit ; ns++ )
1816 *(shifted_slit+ns) = *(slit+ns) ;
1817 for( nl = 0 ; nl<nlambda ; nl++ )
1818 *(shifted_lambda+nl) = *(lambda+nl) ;
1819
1820 // Calculate nshift from BA_pos, AB_pos and nslit
1821 step_slit = get_step_slit( slit, nslit ) ;
1822
1823 // Now fill the rectified shifted image (data, errs and qual)
1824 // Initialize qual with "QFLAG_OUT_OF_NOD" (TBV)
1825 for( ns = 0 ; ns < nslit ; ns++ )
1826 for( nl = 0 ; nl < nlambda ; nl++ )
1827 *(shifted_qual1+nl+ns*nlambda) = QFLAG_OUT_OF_NOD ;
1828
1829 for( ns = 0 ; ns < nslit ; ns++ ) {
1830 int shifted_idx =0;
1831
1832 for( nl = 0 ; nl<nlambda ; nl++ ) {
1833 double ab_slit, ba_slit ;
1834
1835 ab_slit = cpl_polynomial_eval_1d( loc0_list->cenpoly,
1836 *(lambda+nl), NULL ) ;
1837 ba_slit = cpl_polynomial_eval_1d( loc_list->cenpoly,
1838 *(lambda+nl), NULL ) ;
1839 nshift = (float)(ab_slit-ba_slit)/step_slit ;
1840 xsh_msg_dbg_medium( " nl: %d, ns: %d, nshift: %d, ab: %lf, ba: %lf",
1841 nl, ns, nshift, ab_slit, ba_slit ) ;
1842
1843 shifted_idx = nl+(ns+nshift)*nlambda ;
1844
1845 if ( shifted_idx < 0 || shifted_idx > max_size ) {
1846 xsh_msg_dbg_high( " Out of Bound: nl=%d, ns=%d, shifted=%d, max=%d", nl, ns,
1847 shifted_idx, max_size ) ;
1848 break ;
1849 }
1850 *(shifted_flux1+nl+(ns+nshift)*nlambda) =
1851 *(flux1+nl+ns*nlambda) ;
1852 *(shifted_errs1+nl+(ns+nshift)*nlambda) =
1853 *(errs1+nl+ns*nlambda) ;
1854 *(shifted_qual1+nl+(ns+nshift)*nlambda) =
1855 *(qual1+nl+ns*nlambda) ;
1856 }
1857 if ( shifted_idx > max_size ) break ;
1858 }
1859 }
1860
1861 // Save shifted frame
1862 tag = xsh_stringcat_any( cpl_frame_get_tag( rec_frame ), "", (void*)NULL );
1863 tag_drl = xsh_stringcat_any( cpl_frame_get_tag( rec_frame ), "_DRL",
1864 (void*)NULL ) ;
1865 fname_drl = xsh_stringcat_any( tag_drl, ".fits", (void*)NULL ) ;
1866 check(*res_frame_ext=xsh_rec_list_save2(shifted_rec_list,fname,tag)) ;
1867
1868
1869 check( result =xsh_rec_list_save( shifted_rec_list, fname_drl,
1870 tag_drl, CPL_TRUE));
1871
1872 xsh_add_temporary_file(fname_drl);
1873
1874 xsh_msg( " Shifted Frame SAVED [%s]", cpl_frame_get_tag( result ) ) ;
1875
1876 cleanup:
1877 xsh_localization_free( &loc_list ) ;
1878 xsh_localization_free( &loc0_list ) ;
1879 xsh_rec_list_free( &rec_list ) ;
1880 xsh_rec_list_free( &shifted_rec_list ) ;
1881
1882 return result ;
1883}
1884*/
1885
1886static double
1887compute_shift_with_kw( cpl_propertylist *header,
1888 xsh_rectify_param *rectify_par,
1889 double **ref_ra,
1890 double **ref_dec,
1891 int flag)
1892{
1893 double ref_ra_cumoff, ref_dec_cumoff;
1894 double ra_cumoff, dec_cumoff;
1895 double ra_reloff, dec_reloff, ra_off, dec_off;
1896 double posang=0;
1897 /* double bin_space=0; */
1898 double shift_arcsec=0;
1899
1900
1901
1902 xsh_msg_dbg_high( "==> compute_shift_with_kw" ) ;
1903 XSH_ASSURE_NOT_NULL( header);
1904 XSH_ASSURE_NOT_NULL( rectify_par);
1905
1906 check( posang = xsh_pfits_get_posang( header));
1907 posang = posang/180.0* M_PI;
1908
1909 if ( *ref_ra != NULL){
1910 ref_ra_cumoff = **ref_ra;
1911 }
1912 else{
1913 check( ref_ra_cumoff = xsh_pfits_get_ra_cumoffset( header));
1914 XSH_MALLOC( *ref_ra, double, 1);
1915 **ref_ra = ref_ra_cumoff;
1916 }
1917 if ( *ref_dec != NULL){
1918 ref_dec_cumoff = **ref_dec;
1919 }
1920 else{
1921 check( ref_dec_cumoff = xsh_pfits_get_dec_cumoffset( header));
1922 XSH_MALLOC( *ref_dec, double, 1);
1923 **ref_dec = ref_dec_cumoff;
1924 }
1925
1926 if ( flag == 0){
1927 check( ra_cumoff = xsh_pfits_get_ra_cumoffset( header));
1928 check( dec_cumoff = xsh_pfits_get_dec_cumoffset( header));
1929 }
1930 else{
1931 check( ra_cumoff = xsh_pfits_get_b_ra_cumoffset( header));
1932 check( dec_cumoff = xsh_pfits_get_b_dec_cumoffset( header));
1933 }
1934
1935 check( ra_reloff = xsh_pfits_get_b_ra_reloffset( header));
1936 check( dec_reloff = xsh_pfits_get_b_dec_reloffset( header));
1937
1938 /* compute the move in slit */
1939 xsh_msg( "POSANG %f (rad) REF CUM_(RA DEC) : %f %f ", posang, ref_ra_cumoff,
1940 ref_dec_cumoff);
1941 xsh_msg( "OBJ CUM_(RA DEC) : %f %f ", ra_cumoff, dec_cumoff);
1942 xsh_msg( "REL_(RA DEC) : %f %f", ra_reloff, dec_reloff);
1943
1944 ra_off = ra_cumoff-ref_ra_cumoff;
1945 dec_off = dec_cumoff-ref_dec_cumoff;
1946
1947 xsh_msg( "COMPUTE OFF_(RA DEC) : %f %f", ra_off, dec_off);
1948
1949 shift_arcsec = cos(-posang)*dec_off+sin(-posang)*ra_off;
1950
1951 xsh_msg( "COMPUTE shift in arcsec %f :", shift_arcsec);
1952
1953 cleanup:
1954 return shift_arcsec;
1955}
1956
1957/*****************************************************************************/
1958/*****************************************************************************/
1973cpl_frame*
1974shift_with_kw( cpl_frame *rec_frame,
1976 xsh_rectify_param *rectify_par,
1977 const char *fname,
1978 cpl_frame** res_frame_ext,
1979 double **ref_ra,
1980 double **ref_dec,
1981 int flag)
1982{
1983 cpl_frame *result = NULL ;
1984 int shift_pix;
1985 double bin_space, shift_arcsec;
1986 const char *filename = NULL;
1987 cpl_propertylist* header = NULL;
1988
1989
1990 XSH_ASSURE_NOT_NULL( rec_frame);
1992 XSH_ASSURE_NOT_NULL( fname);
1993 XSH_ASSURE_NOT_NULL( res_frame_ext);
1994 XSH_ASSURE_NOT_NULL( ref_ra);
1995 XSH_ASSURE_NOT_NULL( ref_dec);
1996
1997 check( filename = cpl_frame_get_filename( rec_frame));
1998 check( header = cpl_propertylist_load( filename, 0));
1999
2000 check( bin_space = xsh_pfits_get_rectify_bin_space( header));
2001
2002 check( shift_arcsec = compute_shift_with_kw( header,
2003 rectify_par, ref_ra, ref_dec, flag));
2004
2005 /* here we compute integer shifts not to resample twice in fast mode */
2006 shift_pix = xsh_round_double(shift_arcsec/bin_space);
2007
2008 shift_arcsec = shift_pix*bin_space;
2009
2010 xsh_msg( "SHIFT A-->B : %f in arcsec", shift_arcsec);
2011
2012 check( result = xsh_shift( rec_frame, instrument, fname, shift_arcsec,
2013 res_frame_ext));
2014
2015 cleanup:
2016 if ( cpl_error_get_code() !=CPL_ERROR_NONE){
2017 xsh_free_frame( &result);
2018 }
2019 xsh_free_propertylist( &header);
2020 return result ;
2021}
2022/*****************************************************************************/
2023
2024/*****************************************************************************/
2025/*****************************************************************************/
2026static cpl_frame*
2027xsh_shift( cpl_frame *rec_frame,
2029 const char *fname,
2030 double slit_shift,
2031 cpl_frame** res_frame_ext)
2032{
2033 xsh_rec_list *rec_list = NULL ;
2034 cpl_frame *result = NULL ;
2035 int nb_orders, order = 0 ;
2036 float *slit = NULL ;
2037 //int shift_pix;
2038 char* fname_drl=NULL;
2039 char* tag_drl=NULL;
2040 char* tag=NULL;
2041 int nslit=0;
2042
2043 XSH_ASSURE_NOT_NULL( rec_frame);
2045 XSH_ASSURE_NOT_NULL( fname);
2046 XSH_ASSURE_NOT_NULL( res_frame_ext);
2047
2048 check( rec_list = xsh_rec_list_load( rec_frame, instrument));
2049
2050 nb_orders = rec_list->size ;
2051 check( nslit = xsh_rec_list_get_nslit( rec_list, 0));
2052
2053 /* Loop over orders */
2054 for( order = 0 ; order < nb_orders ; order++ ) {
2055 int ns;
2056
2057 check( slit = xsh_rec_list_get_slit( rec_list, order));
2058 for( ns = 0 ; ns < nslit ; ns++ ) {
2059 slit[ns]+=slit_shift;
2060 }
2061 }
2062 rec_list->slit_min = slit[0];
2063 rec_list->slit_max = slit[nslit-1];
2064
2066 rec_list->slit_min));
2068 rec_list->slit_max));
2069
2070 /* Save shifted frame */
2071 tag= xsh_stringcat_any( cpl_frame_get_tag( rec_frame ), "", (void*)NULL ) ;
2072 tag_drl = xsh_stringcat_any( cpl_frame_get_tag( rec_frame ), "_DRL",
2073 (void*)NULL ) ;
2074 fname_drl = xsh_stringcat_any( "DRL_", fname, (void*)NULL);
2075
2076 check( *res_frame_ext=xsh_rec_list_save2( rec_list, fname,tag)) ;
2077 check( result = xsh_rec_list_save( rec_list, fname_drl,
2078 tag_drl, CPL_TRUE));
2079 xsh_add_temporary_file(fname_drl);
2080
2081 cleanup:
2082 if ( cpl_error_get_code() !=CPL_ERROR_NONE){
2083 xsh_free_frame( &result);
2084 }
2085 XSH_FREE( fname_drl);
2086 XSH_FREE( tag_drl);
2087 XSH_FREE( tag);
2088 xsh_rec_list_free( &rec_list);
2089 return result ;
2090}
2091
2108cpl_frame * xsh_shift_rectified( cpl_frame * rec_frame,
2109 cpl_frame * loc_frame,
2110 cpl_frame * loc0_frame,
2111 const char * file_name,
2112 xsh_combine_nod_param * combine_nod_param,
2113 xsh_rectify_param * rectif_par,
2115 cpl_frame** res_frame_ext )
2116{
2117 cpl_frame * result = NULL ;
2118 double *ra = NULL, *dec = NULL;
2119
2120 xsh_msg( "===> xsh_shift_rectified" ) ;
2121
2122 XSH_ASSURE_NOT_NULL( combine_nod_param ) ;
2123 XSH_ASSURE_NOT_NULL( rec_frame ) ;
2124
2125 check( result = shift_with_kw( rec_frame, instrument,
2126 rectif_par, file_name, res_frame_ext, &ra, &dec, 1));
2127
2128 cleanup:
2129 return result ;
2130}
2131/*----------------------------------------------------------------------------*/
2132
static const double step
static char mode[32]
static double sigma
static xsh_instrument * instrument
static float slit_step
cpl_mask * xsh_qual_to_cpl_mask(cpl_image *qual, const int decode_bp)
void xsh_baryvel(const cpl_propertylist *raw_header, double *bary_corr, double *helio_corr)
Compute velocity correction.
Definition: xsh_baryvel.c:95
xsh_dispersol_list * xsh_dispersol_list_load(cpl_frame *frame, xsh_instrument *instrument)
Load a dispersion list from a frame.
void xsh_dispersol_list_free(xsh_dispersol_list **list)
Free the dispersion list.
xsh_localization * xsh_localization_load(cpl_frame *frame)
Load a localization list from a frame.
void xsh_localization_free(xsh_localization **list)
free memory associated to a localization_list
xsh_order_list * xsh_order_list_load(cpl_frame *frame, xsh_instrument *instr)
load an order list from a frame
void xsh_order_list_free(xsh_order_list **list)
free memory associated to an order_list
cpl_image * xsh_pre_get_data(xsh_pre *pre)
Get data.
xsh_pre * xsh_pre_load(cpl_frame *frame, xsh_instrument *instr)
Load a xsh_pre structure from a frame.
Definition: xsh_data_pre.c:849
cpl_image * xsh_pre_get_qual(xsh_pre *pre)
Get qual.
int xsh_pre_get_ny(const xsh_pre *pre)
Get ny of pre structure.
void xsh_pre_free(xsh_pre **pre)
Free a xsh_pre structure.
Definition: xsh_data_pre.c:823
int xsh_pre_get_nx(const xsh_pre *pre)
Get nx of pre structure.
cpl_image * xsh_pre_get_errs(xsh_pre *pre)
Get errs.
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)
cpl_frame * xsh_rec_list_save(xsh_rec_list *list, const char *filename, const char *tag, int is_temp)
Save a rec list in a frame.
cpl_frame * xsh_rec_list_save2(xsh_rec_list *list, const char *filename, const char *tag)
save an rec list to a frame
float * xsh_rec_list_get_errs1(xsh_rec_list *list, int idx)
void xsh_rec_list_update_header(xsh_rec_list *rec_list, xsh_pre *pre, xsh_rectify_param *rec_par, const char *pro_catg)
Update header of rectified list writing mandatory KW.
int xsh_rec_list_get_order(xsh_rec_list *list, int idx)
Definition: xsh_data_rec.c:846
void xsh_rec_list_set_data_size(xsh_rec_list *list, int idx, int absorder, int nlambda, int ns)
Allocate memory for the order idx of the rectify list.
Definition: xsh_data_rec.c:191
int * xsh_rec_list_get_qual1(xsh_rec_list *list, int idx)
float * xsh_rec_list_get_slit(xsh_rec_list *list, int idx)
Definition: xsh_data_rec.c:884
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
xsh_rec_list * xsh_rec_list_create(xsh_instrument *instr)
Create an empty order list.
Definition: xsh_data_rec.c:100
cpl_frame * xsh_rec_list_save_table(xsh_rec_list *list, const char *filename, const char *tag, int is_temp)
Save a rec list in a frame.
xsh_spectralformat_list * xsh_spectralformat_list_load(cpl_frame *frame, xsh_instrument *instr)
Load a spectralformat list from a frame.
void xsh_spectralformat_list_free(xsh_spectralformat_list **list)
Free memory associated to an spactralformat_list.
xsh_wavesol * xsh_wavesol_load(cpl_frame *frame, xsh_instrument *instrument)
load a wavelength solution
double xsh_wavesol_eval_poly(xsh_wavesol *sol, double lambda, double order, double slit)
eval the polynomial solution in Y
void xsh_wavesol_free(xsh_wavesol **w)
free wavelength solution structure
double xsh_wavesol_eval_polx(xsh_wavesol *sol, double lambda, double order, double slit)
eval the polynomial solution in X
#define XSH_REGDEBUG(...)
Definition: xsh_error.h:132
#define XSH_ASSURE_NOT_ILLEGAL(cond)
Definition: xsh_error.h:107
#define assure(CONDITION, ERROR_CODE,...)
Definition: xsh_error.h:54
#define check(COMMAND)
Definition: xsh_error.h:71
#define XSH_CMP_INT(A, OPERATOR, B, SUFFIX,...)
Definition: xsh_error.h:119
#define xsh_error_reset()
Definition: xsh_error.h:87
#define xsh_error_msg(...)
Definition: xsh_error.h:94
#define XSH_ASSURE_NOT_NULL(pointer)
Definition: xsh_error.h:99
const char * xsh_instrument_arm_tostring(xsh_instrument *i)
Get the string associated with an arm.
XSH_MODE xsh_instrument_get_mode(xsh_instrument *i)
Get a mode on instrument structure.
cpl_error_code xsh_model_config_load_best(cpl_frame *config_frame, xsh_xs_3 *p_xs_3)
Load the config model table and fill the struct.
Definition: xsh_model_io.c:174
void xsh_model_get_xy(xsh_xs_3 *p_xs_3, xsh_instrument *instr, double lambda_nm, int morder, double ent_slit_pos, double *x, double *y)
Compute the detector location (floating point pixels) of a given wavelength/entrance slit position.
int size
int * y
#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(...)
Print a message on info level.
Definition: xsh_msg.h:121
#define xsh_msg_dbg_low(...)
Definition: xsh_msg.h:48
#define xsh_msg_dbg_high(...)
Definition: xsh_msg.h:40
double xsh_pfits_get_ra_cumoffset(const cpl_propertylist *plist)
Definition: xsh_pfits.c:3585
double xsh_pfits_get_b_ra_reloffset(const cpl_propertylist *plist)
Definition: xsh_pfits.c:3615
double xsh_pfits_get_shiftifu_lambdaref(cpl_propertylist *plist)
Definition: xsh_pfits.c:4030
double xsh_pfits_get_slitmap_median_edglo(const cpl_propertylist *plist)
Definition: xsh_pfits.c:3914
void xsh_pfits_set_shiftifu_lambdaref(cpl_propertylist *plist, double value)
Definition: xsh_pfits.c:4015
void xsh_pfits_set_rectify_space_max(cpl_propertylist *plist, double value)
WRITE the space (slit) max value.
Definition: xsh_pfits.c:3268
double xsh_pfits_get_rectify_bin_space(cpl_propertylist *plist)
find out the rectify space (slit) binning
Definition: xsh_pfits.c:3302
double xsh_pfits_get_b_ra_cumoffset(const cpl_propertylist *plist)
Definition: xsh_pfits.c:3645
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_rectify_space_min(cpl_propertylist *plist, double value)
WRITE the space (slit) min value.
Definition: xsh_pfits.c:3252
double xsh_pfits_get_posang(const cpl_propertylist *plist)
Definition: xsh_pfits.c:3513
double xsh_pfits_get_dec_cumoffset(const cpl_propertylist *plist)
Definition: xsh_pfits.c:3600
double xsh_pfits_get_slitmap_median_slicup(const cpl_propertylist *plist)
Definition: xsh_pfits.c:3946
double xsh_pfits_get_slitmap_median_edgup(const cpl_propertylist *plist)
Definition: xsh_pfits.c:3898
double xsh_pfits_get_slitmap_median_sliclo(const cpl_propertylist *plist)
Definition: xsh_pfits.c:3962
double xsh_pfits_get_b_dec_cumoffset(const cpl_propertylist *plist)
Definition: xsh_pfits.c:3660
double xsh_pfits_get_b_dec_reloffset(const cpl_propertylist *plist)
Definition: xsh_pfits.c:3630
static double compute_shift_with_localization(cpl_frame *loc_frame, cpl_frame *loc0_frame)
Definition: xsh_rectify.c:1708
cpl_frame * xsh_rectify_and_shift(cpl_frame *sci_frame, cpl_frame *orderlist_frame, cpl_frame *wavesol_frame, cpl_frame *model_frame, xsh_instrument *instrument, xsh_rectify_param *rectify_par, cpl_frame *spectralformat_frame, cpl_frame *loc_frame, cpl_frame *loc0_frame, double *throw_shift, cpl_frame *disp_tab_frame, const char *res_name, cpl_frame **res_frame_ext, cpl_frame **res_frame_tab)
Definition: xsh_rectify.c:928
cpl_frameset * xsh_rectify_orders_ifu(cpl_frame *sci_frame, xsh_order_list *orderlist, cpl_frameset *wavesol_frameset, cpl_frameset *shift_frameset, cpl_frame *model_frame, xsh_instrument *instrument, xsh_rectify_param *rectify_par, cpl_frame *spectralformat_frame, cpl_frame *slitmap_frame, cpl_frameset **res_frameset_ext, cpl_frameset **res_frameset_tab, int min_index, int max_index, const char *rec_prefix)
Definition: xsh_rectify.c:1514
cpl_frame * shift_with_kw(cpl_frame *rec_frame, xsh_instrument *instrument, xsh_rectify_param *rectify_par, const char *fname, cpl_frame **res_frame_ext, double **ref_ra, double **ref_dec, int flag)
This function creates a structure containing for each order the shift to be applied.
Definition: xsh_rectify.c:1974
cpl_frame * xsh_rectify_orders(cpl_frame *sci_frame, xsh_order_list *orderlist, cpl_frame *wavesol_frame, cpl_frame *model_frame, xsh_instrument *instrument, xsh_rectify_param *rectify_par, cpl_frame *spectralformat_frame, cpl_frame *disp_tab_frame, const char *res_name, const char *tag, cpl_frame **res_frame_ext, cpl_frame **res_frame_tab, int min_index, int max_index, double slit_min, int nslit, double slit_shift, cpl_frame *slitshift_tab)
Create a grid in wavelength with the step in lambda and slit. Steps are defined in the parameters....
Definition: xsh_rectify.c:1079
void xsh_compute_slitlet_limits(cpl_frameset *shift_set, double sdown, double sldown, double slup, double sup, double slit_bin, double *slitmin_tab, int *nslit_tab, double *slitcen_tab)
Definition: xsh_rectify.c:1372
void xsh_get_slit_edges(cpl_frame *slitmap_frame, double *sdown, double *sup, double *sldown, double *slup, xsh_instrument *instrument)
Trace slit edges in a master flat.
Definition: xsh_rectify.c:671
void xsh_rec_slit_size(xsh_rectify_param *rectify_par, double *slit_min, int *nslit, XSH_MODE mode)
rectify frame
Definition: xsh_rectify.c:758
cpl_frame * xsh_rectify(cpl_frame *sci_frame, cpl_frame *orderlist_frame, cpl_frame *wavesol_frame, cpl_frame *model_frame, xsh_instrument *instrument, xsh_rectify_param *rectify_par, cpl_frame *spectralformat_frame, cpl_frame *disp_tab_frame, const char *res_name, cpl_frame **res_frame_ext, cpl_frame **res_frame_tab, const char *rec_prefix)
Create a grid in wavelength with the step in lambda and slit. Steps are defined in the parameters....
Definition: xsh_rectify.c:867
static void xsh_rec_list_rectify(xsh_rec_list *rec_list, int iorder, int irec, xsh_pre *sci_pre, cpl_image *var_img, double fx, double fy, double radius, cpl_vector *profile, cpl_vector *profile2, double mult)
Definition: xsh_rectify.c:188
static void adjust_lambdas(xsh_spectralformat_list *spec_list, xsh_rectify_param *rectify_par)
Definition: xsh_rectify.c:788
static cpl_vector * rec_profile
Definition: xsh_rectify.c:80
static cpl_frame * xsh_shift(cpl_frame *rec_frame, xsh_instrument *instrument, const char *fname, double slit_shift, cpl_frame **res_frame_ext)
Definition: xsh_rectify.c:2027
static void fill_rectified(xsh_pre *pre_sci, xsh_rec_list *rec_list, int idx, xsh_wavesol *wavesol, xsh_xs_3 *model_config, xsh_instrument *instrument, xsh_dispersol_list *disp_list, float slit_min, float slit_max, double lambda_min, int skip_low, int skip_up, xsh_rectify_param *rectify_par, double slit_shift, cpl_frame *slit_shiftab_frame)
Fill the rectified structure for one order with the corresponding flux, using the wavelength solution...
Definition: xsh_rectify.c:371
static void xsh_frame_set_shiftifu_ref(cpl_propertylist *header, cpl_frame *shift_frame)
Definition: xsh_rectify.c:1332
cpl_frameset * xsh_rectify_ifu(cpl_frame *sci_frame, cpl_frame *orderlist_frame, cpl_frameset *wavesol_frameset, cpl_frameset *shiftifu_frameset, cpl_frame *model_frame, xsh_instrument *instrument, xsh_rectify_param *rectify_par, cpl_frame *spectralformat_frame, cpl_frame *slitmap_frame, cpl_frameset **rec_frameset_ext, cpl_frameset **rec_frameset_tab, const char *rec_prefix)
Definition: xsh_rectify.c:1014
cpl_frame * xsh_shift_rectified(cpl_frame *rec_frame, cpl_frame *loc_frame, cpl_frame *loc0_frame, const char *file_name, xsh_combine_nod_param *combine_nod_param, xsh_rectify_param *rectif_par, xsh_instrument *instrument, cpl_frame **res_frame_ext)
Definition: xsh_rectify.c:2108
static double compute_shift_with_kw(cpl_propertylist *header, xsh_rectify_param *rectify_par, double **ref_ra, double **ref_dec, int flag)
Definition: xsh_rectify.c:1887
static cpl_vector * err_profile
Definition: xsh_rectify.c:81
static void xsh_get_shift_ref(cpl_frameset *set, double *down, double *cen, double *up)
Definition: xsh_rectify.c:1291
void xsh_free_polynomial(cpl_polynomial **p)
Deallocate a polynomial and set the pointer to NULL.
Definition: xsh_utils.c:2194
void xsh_free_vector(cpl_vector **v)
Deallocate a vector and set the pointer to NULL.
Definition: xsh_utils.c:2284
void xsh_free_image(cpl_image **i)
Deallocate an image and set the pointer to NULL.
Definition: xsh_utils.c:2116
int xsh_debug_level_get(void)
get debug level
Definition: xsh_utils.c:3142
double convert_data_to_bin(double data, int binning)
Definition: xsh_utils.c:3246
void xsh_free_frame(cpl_frame **f)
Deallocate a frame and set the pointer to NULL.
Definition: xsh_utils.c:2269
double xsh_data_interpolate(double wav, int nrow, double *pw, double *pe)
Interpolate data points.
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
void xsh_add_temporary_file(const char *name)
Add temporary file to temprary files list.
Definition: xsh_utils.c:1432
xsh_dispersol * list
cpl_polynomial * slit_poly
cpl_polynomial * lambda_poly
cpl_polynomial * cenpoly
cpl_image * qual
Definition: xsh_data_pre.h:71
cpl_image * data
Definition: xsh_data_pre.h:65
cpl_propertylist * data_header
Definition: xsh_data_pre.h:66
double slit_max
Definition: xsh_data_rec.h:83
cpl_propertylist * header
Definition: xsh_data_rec.h:89
double slit_min
Definition: xsh_data_rec.h:82
xsh_rec * list
Definition: xsh_data_rec.h:87
float * slit
Definition: xsh_data_rec.h:66
int nslit
Definition: xsh_data_rec.h:65
int nlambda
Definition: xsh_data_rec.h:64
double * lambda
Definition: xsh_data_rec.h:67
cpl_kernel kernel_type
#define QFLAG_OUTSIDE_DATA_RANGE
#define QFLAG_MISSING_DATA
@ XSH_MODE_SLIT
@ XSH_MODE_IFU
#define MIN_SLIT
#define MAX_SLIT
int nx
int n
Definition: xsh_detmon_lg.c:92
int ny
int order
Definition: xsh_detmon_lg.c:80
#define XSH_ORDER2D
Definition: xsh_dfs.h:588
#define XSH_GET_TAG_FROM_ARM(TAG, instr)
Definition: xsh_dfs.h:1548
#define XSH_RECTIFY_TYPE_MODEL
Definition: xsh_drl.h:378
#define XSH_SHIFTIFU_COLNAME_WAVELENGTH
Definition: xsh_drl.h:516
#define XSH_SHIFTIFU_COLNAME_SHIFTSLIT
Definition: xsh_drl.h:518
#define XSH_RECTIFY_TYPE_POLY
Definition: xsh_drl.h:377
#define max(a, b)
@ CENTER_IFU_SLITLET
Definition: xsh_ifu_defs.h:42
@ LOWER_IFU_SLITLET
Definition: xsh_ifu_defs.h:42
@ UPPER_IFU_SLITLET
Definition: xsh_ifu_defs.h:42
#define SLITLET_CEN_CENTER
Definition: xsh_ifu_defs.h:36
cpl_error_code xsh_model_temperature_update_frame(cpl_frame **model_config_frame, cpl_frame *ref_frame, xsh_instrument *instrument, int *found_temp)
#define XSH_QC_VRAD_HELICOR_C
Definition: xsh_pfits_qc.h:47
#define XSH_QC_VRAD_HELICOR
Definition: xsh_pfits_qc.h:46
#define XSH_QC_VRAD_BARYCOR
Definition: xsh_pfits_qc.h:44
#define XSH_QC_VRAD_BARYCOR_C
Definition: xsh_pfits_qc.h:45
#define XSH_FREE(POINTER)
Definition: xsh_utils.h:92
@ XSH_DEBUG_LEVEL_MEDIUM
Definition: xsh_utils.h:138
#define XSH_MALLOC(POINTER, TYPE, SIZE)
Definition: xsh_utils.h:49
#define XSH_CALLOC(POINTER, TYPE, SIZE)
Definition: xsh_utils.h:56
#define M_PI
Definition: xsh_utils.h:43
#define XSH_TABLE_LOAD(TABLE, NAME)
#define XSH_TABLE_FREE(TABLE)