X-shooter Pipeline Reference Manual 3.8.15
xsh_flexcor.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-10-02 11:35:51 $
23 * $Revision: 1.30 $
24 * $Name: not supported by cvs2svn $
25 */
26
27#ifdef HAVE_CONFIG_H
28#include <config.h>
29#endif
30
31/*---------------------------------------------------------------------------*/
37/*---------------------------------------------------------------------------*/
40/*----------------------------------------------------------------------------
41 Includes
42 ----------------------------------------------------------------------------*/
43
44#include <xsh_drl.h>
45#include <xsh_pfits.h>
46#include <xsh_data_resid_tab.h>
47#include <xsh_data_instrument.h>
48#include <xsh_data_arclist.h>
50#include <xsh_data_order.h>
51#include <xsh_data_wavesol.h>
52#include <xsh_data_the_map.h>
53#include <xsh_model_io.h>
54#include <xsh_model_utils.h>
55#include <xsh_error.h>
56#include <xsh_msg.h>
57#include <cpl.h>
58
59/*----------------------------------------------------------------------------
60 Typedefs
61 ---------------------------------------------------------------------------*/
62
63/*----------------------------------------------------------------------------
64 Functions prototypes
65 ---------------------------------------------------------------------------*/
66
67
68/*----------------------------------------------------------------------------
69 Implementation
70 ---------------------------------------------------------------------------*/
71
72
73/*---------------------------------------------------------------------------*/
90/*---------------------------------------------------------------------------*/
91cpl_frame* xsh_flexcor( cpl_frame* afc_frame,
92 cpl_frame *wave_tab_frame,
93 cpl_frame *model_config_frame,
94 cpl_frame *order_tab_frame,
95 cpl_frame *attresidtab_frame,
96 int afc_xmin, int afc_ymin,
97 xsh_instrument *instr,
98 cpl_frame **afc_order_tab_frame,
99 cpl_frame **afc_model_config_frame)
100{
101 char result_name[256];
102 cpl_frame *result = NULL;
103 cpl_frame *afc_order_tab_result = NULL;
104 cpl_table *trace = NULL;
105 float xshift = 0.0, yshift = 0.0;
106 float T_effect_x = 0.0, T_effect_y = 0.0;
107 xsh_resid_tab *attresid_tab = NULL;
108 int residtab_size, iline;
109 const char *tag = NULL;
110 xsh_wavesol *wavesol = NULL;
111 xsh_order_list *order_tab = NULL;
112 double lambda, order, slit=0.0;
113 double att_x, att_y;
114 double cal_x, cal_y;
115 double cal_no_T_corr_x, cal_no_T_corr_y;
116 XSH_INSTRCONFIG *config = NULL;
118 xsh_xs_3 config_model;
119 xsh_xs_3 cfg_mod_no_T_corr;
120 // xsh_xs_3* p_xs_3;
121 cpl_table *model_tab = NULL;
122 cpl_propertylist *model_header = NULL;
123 cpl_propertylist *afc_header = NULL;
124 cpl_vector *shiftx_vect = NULL;
125 cpl_vector *shifty_vect = NULL;
126 cpl_vector *T_effect_x_vect = NULL;
127 cpl_vector *T_effect_y_vect = NULL;
128 const char* afc_name=NULL;
129 int found_temp=true;
130 cpl_frame* MODEL_CONF_OPT_frame=NULL;
131 cpl_propertylist* qc_head2copy=NULL;
132 //xsh_pre * pre_afc = NULL ;
133 //xsh_msg("p1=%p p2=%p",afc_order_tab_frame,*afc_model_config_frame);
134 /* check input parameters */
135 XSH_ASSURE_NOT_NULL( attresidtab_frame);
136 XSH_ASSURE_NOT_NULL( order_tab_frame);
137 XSH_ASSURE_NOT_NULL( afc_order_tab_frame);
138 //XSH_ASSURE_NOT_NULL( afc_model_config_frame);
139
140 /* First compute the shift */
141 check( attresid_tab = xsh_resid_tab_load( attresidtab_frame));
142 check( order_tab = xsh_order_list_load( order_tab_frame, instr));
143 check( residtab_size = xsh_resid_tab_get_size( attresid_tab));
144
145 XSH_REGDEBUG( "RESID_TAB SIZE %d", residtab_size);
146 check(afc_name=cpl_frame_get_filename(afc_frame));
147 check(afc_header=cpl_propertylist_load(afc_name,0));
148
149 if ( wave_tab_frame != NULL){
150 check( wavesol = xsh_wavesol_load( wave_tab_frame, instr));
151 }
152 else{
153 XSH_ASSURE_NOT_NULL( model_config_frame);
154 check( xsh_model_config_load_best( model_config_frame, &cfg_mod_no_T_corr));
155 check(xsh_model_temperature_update_frame(&model_config_frame,afc_frame,
156 instr,&found_temp));
157 check( xsh_model_config_load_best( model_config_frame, &config_model));
158 }
159
160 check( shiftx_vect = cpl_vector_new( residtab_size));
161 check( shifty_vect = cpl_vector_new( residtab_size));
162 check( T_effect_x_vect = cpl_vector_new( residtab_size));
163 check( T_effect_y_vect = cpl_vector_new( residtab_size));
164
165 check(arm=xsh_instrument_get_arm(instr));
166
167 for( iline=0; iline < residtab_size; iline++){
168
169 lambda = attresid_tab->lambda[iline];
170 order = attresid_tab->order[iline];
171 att_x = attresid_tab->xgauss[iline];
172 att_y = attresid_tab->ygauss[iline];
173
174 if ( wavesol != NULL){
175 check( cal_x = xsh_wavesol_eval_polx( wavesol, lambda, order, slit));
176 check( cal_y = xsh_wavesol_eval_poly( wavesol, lambda, order, slit));
177 }
178 else{
179 check( xsh_model_get_xy( &config_model, instr, lambda, order, slit,
180 &cal_x, &cal_y));
181 check( xsh_model_get_xy( &cfg_mod_no_T_corr, instr, lambda, order, slit,
182 &cal_no_T_corr_x, &cal_no_T_corr_y));
183 /*thpre_x/y are only necessary for the pre-annealing residuals check
184 in xsh_model_pipe_anneal*/
185 attresid_tab->thpre_x[iline]=cal_x;
186 attresid_tab->thpre_y[iline]=cal_y;
187
188 attresid_tab->slit_index[iline]=4;
189
190 T_effect_x=cal_x-cal_no_T_corr_x;
191 T_effect_y=cal_y-cal_no_T_corr_y;
192 xsh_msg(" T_effect_x %f T_effect_y %f", T_effect_x,T_effect_y);
193 check( cpl_vector_set( T_effect_x_vect, iline, T_effect_x));
194 check( cpl_vector_set( T_effect_y_vect, iline, T_effect_y));
195 }
196
197 cal_x = cal_x-afc_xmin+1.0;
198 cal_y =cal_y-afc_ymin+1.0;
199
200 xshift = att_x-cal_x;
201 yshift = att_y-cal_y;
202 xsh_msg("lambda %f order %f slit %f calx %f caly %f shiftx %f shifty %f TeffX %f TeffY %f",
203 lambda,order,slit,cal_x, cal_y,xshift,yshift, T_effect_x,T_effect_y);
204
205 check( cpl_vector_set( shiftx_vect, iline, xshift));
206 check( cpl_vector_set( shifty_vect, iline, yshift));
207 }
208
209 check( xshift = cpl_vector_get_median( shiftx_vect));
210 check( yshift = cpl_vector_get_median( shifty_vect));
211 check( T_effect_x = cpl_vector_get_median( T_effect_x_vect));
212 check( T_effect_y = cpl_vector_get_median( T_effect_y_vect));
213
214 xsh_msg("Measured shift IN X %f Y %f", xshift, yshift);
215 if ( wavesol == NULL){
216 xsh_msg("Temperature effect at (mean) reference line in X %f Y %f",
217 T_effect_x,T_effect_y);
218 }
219
220
221 if ( wavesol != NULL){
222 cpl_table* tab=NULL;
223 cpl_propertylist* plist=NULL;
224 check( xsh_wavesol_apply_shift( wavesol, xshift, yshift));
225 check( trace = xsh_wavesol_trace( wavesol, &lambda, &order, &slit, 1));
227 sprintf( result_name, "%s.fits", tag);
228 check( result = xsh_wavesol_save( wavesol, trace, result_name, tag));
229 tab=cpl_table_load(result_name,1,0);
230 plist=cpl_propertylist_load(result_name,0);
231 cpl_propertylist_append_double(plist,XSH_QC_AFC_XSHIFT,xshift);
232 cpl_propertylist_append_double(plist,XSH_QC_AFC_YSHIFT,yshift);
233 cpl_propertylist_set_comment(plist,XSH_QC_AFC_XSHIFT,XSH_QC_AFC_XSHIFT_C);
234 cpl_propertylist_set_comment(plist,XSH_QC_AFC_YSHIFT,XSH_QC_AFC_YSHIFT_C);
235 check( cpl_table_save(tab,plist, NULL, result_name,CPL_IO_DEFAULT));
236 xsh_free_table(&tab);
237 xsh_free_propertylist(&plist);
238 }
239 else{
240 if (arm!=XSH_ARM_NIR) {
241 xsh_model_offset( xshift, yshift, &config_model);
242 /*New cfg is now appropriate to the prism temperature(s) in the AFC
243 data, when it is written to the table here it will have that new
244 temperature*/
245 }
246 else {
247 /*For NIR arm there are many more lines, so instead of just computing
248 a shift we can anneal a bunch of parameters.*/
249
250 check(MODEL_CONF_OPT_frame=xsh_model_pipe_anneal( model_config_frame,
251 attresidtab_frame,
252 2000,
253 1.0,
254 9,
255 5));
256
257 /* The result of pipe_anneal is a frame, so we need to extract the
258 config_model data structure to be consistent with what was done
259 for the other arms */
260 check( xsh_model_config_load_best(MODEL_CONF_OPT_frame, &config_model));
261
262 }
263
264 xsh_msg("file name=%s",cpl_frame_get_filename(model_config_frame));
265 qc_head2copy=cpl_propertylist_load_regexp(cpl_frame_get_filename(model_config_frame),0,"^ESO QC MODEL",0);
266
267 check( model_tab = xsh_model_io_output_cfg( &config_model));
269 sprintf( result_name, "%s.fits", tag);
270 check( model_header = cpl_propertylist_new());
271 check( xsh_pfits_set_pcatg( model_header, tag));
272 cpl_propertylist_append_double(model_header,XSH_QC_AFC_XSHIFT,xshift);
273 cpl_propertylist_append_double(model_header,XSH_QC_AFC_YSHIFT,yshift);
274 cpl_propertylist_set_comment(model_header,XSH_QC_AFC_XSHIFT,XSH_QC_AFC_XSHIFT_C);
275 cpl_propertylist_set_comment(model_header,XSH_QC_AFC_YSHIFT,XSH_QC_AFC_YSHIFT_C);
276 /*Add extra QC keywords for xsh_flexcomp recipe */
277 if(cpl_propertylist_has(afc_header,"ESO ADA ABSROT START") &&
278 cpl_propertylist_has(afc_header,"ESO ADA ABSROT END")){
279 double absrot_avg = 0.5*(cpl_propertylist_get_double(afc_header,"ESO ADA ABSROT START") +
280 cpl_propertylist_get_double(afc_header,"ESO ADA ABSROT END"));
281 cpl_propertylist_append_double(model_header,"ESO QC ABSROT",absrot_avg);
282 }
283 if(model_tab){
284 for(cpl_size row=0;row<cpl_table_get_nrow(model_tab);row++){
285 const char* fld = cpl_sprintf("%s",cpl_table_get_string(model_tab,XSH_COL_MODEL_CONF_NAME,row));
286 double val = cpl_table_get_double(model_tab,XSH_COL_MODEL_CONF_BEST,row,NULL);
287 if(fld){
288 if(!strcmp(fld,"chipx")){
289 cpl_propertylist_update_double(model_header,"ESO QC CHIP X",val);
290 }
291 if(!strcmp(fld,"chipy")){
292 cpl_propertylist_update_double(model_header,"ESO QC CHIP Y",val);
293 }
294 }
295 }
296 }
297
298 cpl_propertylist_append(model_header,qc_head2copy);
299
300 xsh_free_propertylist(&qc_head2copy);
301 check( cpl_table_save( model_tab, model_header, NULL, result_name,
302 CPL_IO_DEFAULT));
303 check( *afc_model_config_frame = cpl_frame_new());
304 check( cpl_frame_set_filename( *afc_model_config_frame, result_name));
305 check( cpl_frame_set_tag( *afc_model_config_frame, tag));
306 }
307
308 check (xsh_order_list_apply_shift( order_tab, xshift, yshift));
310 sprintf( result_name, "%s.fits", tag);
311
312 check( config = xsh_instrument_get_config( instr));
313 check( afc_order_tab_result = xsh_order_list_save( order_tab, instr,
314 result_name, tag, config->ny));
315
316 *afc_order_tab_frame = afc_order_tab_result;
317
318 cleanup:
319 if( cpl_error_get_code() != CPL_ERROR_NONE) {
320 xsh_free_frame( &result);
321 xsh_free_frame( &afc_order_tab_result);
322
323 }
324 xsh_resid_tab_free( &attresid_tab);
325 xsh_order_list_free( &order_tab);
326 xsh_free_vector(&T_effect_x_vect);
327 xsh_free_vector(&T_effect_y_vect);
328 xsh_free_vector( &shiftx_vect);
329 xsh_free_vector( &shifty_vect);
330 xsh_free_table( &trace);
331 xsh_free_table( &model_tab);
332 xsh_free_propertylist( &model_header);
333 xsh_free_propertylist( &afc_header);
334 if (wavesol==NULL) {
335 if (arm==XSH_ARM_NIR) {
336 xsh_free_frame(&MODEL_CONF_OPT_frame);
337 }
338 }
339 xsh_wavesol_free( &wavesol);
340 //xsh_pre_free( &pre_afc);
341 return result;
342}
343/*---------------------------------------------------------------------------*/
344
345/*---------------------------------------------------------------------------*/
361/*---------------------------------------------------------------------------*/
362cpl_frame*
363xsh_afcthetab_create( cpl_frame *wave_tab_frame,
364 cpl_frame *model_config_frame,
365 int order,
366 cpl_frame *spectralformat_frame,
367 cpl_frame *arclines_frame,
368 int xmin,
369 int ymin,
370 xsh_instrument *instr,
371 const int clean_tmp)
372{
373
374 cpl_frame *result = NULL;
375 xsh_arclist *arclist = NULL;
376 xsh_wavesol *wavesol = NULL;
377 float lambda, slit_position;
378 double x,y;
379 int slit_index;
380 int size, iline, iorder, ithe=0;
381 xsh_the_map* themap = NULL;
382 int the_size=0;
383 xsh_xs_3 config_model;
384 cpl_vector* orders = NULL;
385 int orders_size;
386 xsh_spectralformat_list *spectrallist = NULL;
387 char thename[256];
388
389 /* check input parameters */
390 XSH_ASSURE_NOT_NULL( arclines_frame);
391 XSH_ASSURE_NOT_NULL( instr);
392
393 check( arclist = xsh_arclist_load( arclines_frame));
394
395 check( size = xsh_arclist_get_size( arclist));
396
397 xsh_msg("size %d", size);
398
399 if ( spectralformat_frame != NULL){
400 check( spectrallist = xsh_spectralformat_list_load(
401 spectralformat_frame, instr));
402
403 for( iline=0; iline < size; iline++){
404 check( lambda = xsh_arclist_get_wavelength(arclist, iline));
405 check( orders = xsh_spectralformat_list_get_orders(spectrallist,
406 lambda));
407 if (orders != NULL){
408 check( orders_size = cpl_vector_get_size( orders));
409 the_size += orders_size;
410 }
411 xsh_free_vector( &orders);
412 }
413 }
414 else{
415 the_size = 1;
416 }
417
418
419 xsh_msg("THE SIZE %d", the_size);
420 check( themap = xsh_the_map_create( the_size));
421
422 if (wave_tab_frame != NULL){
423 check( wavesol = xsh_wavesol_load( wave_tab_frame, instr));
424 }
425 else{
426 XSH_ASSURE_NOT_NULL( model_config_frame);
427 check( xsh_model_config_load_best( model_config_frame,
428 &config_model));
429 }
430
431 for( iline=0; iline < size; iline++){
432 check( lambda = xsh_arclist_get_wavelength(arclist, iline));
433 slit_index = 4;
434 slit_position = 0;
435 if ( spectrallist != NULL){
436 check( orders = xsh_spectralformat_list_get_orders( spectrallist,
437 lambda));
438 }
439 else{
440 orders = cpl_vector_new(1);
441 cpl_vector_set( orders, 0, order);
442 }
443 if (orders != NULL){
444 check( orders_size = cpl_vector_get_size( orders));
445
446 for( iorder=0; iorder < orders_size; iorder++){
447 int cur_order;
448
449 check( cur_order = cpl_vector_get( orders, iorder));
450
451 if (wavesol != NULL){
452 check( x = xsh_wavesol_eval_polx( wavesol, lambda, cur_order,
453 slit_position));
454 check( y = xsh_wavesol_eval_poly( wavesol, lambda, cur_order,
455 slit_position));
456 }
457 else{
458 check( xsh_model_get_xy( &config_model, instr, lambda, cur_order,
459 slit_position, &x, &y));
460 }
461 x = x-xmin+1.0;
462 y = y-ymin+1.0;
463 XSH_REGDEBUG("lambda %f order %d %f %f",lambda, cur_order, x, y);
464 check( xsh_the_map_set_arcline( themap, ithe, lambda, cur_order,
465 slit_index, slit_position, x, y));
466 ithe++;
467 }
468 }
469 xsh_free_vector( &orders);
470 }
471 sprintf( thename, "AFC_THEOTAB_%s.fits",
473 check( result = xsh_the_map_save( themap, thename));
474 if(clean_tmp) {
475 xsh_add_temporary_file(thename);
476 }
477
478 cleanup:
479 if( cpl_error_get_code() != CPL_ERROR_NONE) {
480 xsh_free_vector( &orders);
481 xsh_free_frame( &result);
482 }
483 xsh_spectralformat_list_free( &spectrallist);
484 xsh_arclist_free( &arclist);
485 xsh_wavesol_free( &wavesol);
486 xsh_the_map_free( &themap);
487
488 return result;
489}
490/*---------------------------------------------------------------------------*/
float xsh_arclist_get_wavelength(xsh_arclist *list, int idx)
get wavelength of a line in the arcline list
void xsh_arclist_free(xsh_arclist **list)
free memory associated to a arclist
int xsh_arclist_get_size(xsh_arclist *list)
get size of arcline list
xsh_arclist * xsh_arclist_load(cpl_frame *frame)
load an arcline list frame in arclist structure
void xsh_order_list_apply_shift(xsh_order_list *list, double xshift, double yshift)
Shift a order list.
cpl_frame * xsh_order_list_save(xsh_order_list *order_list, xsh_instrument *instrument, const char *filename, const char *tag, const int ny)
Save an order list to a frame.
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
void xsh_resid_tab_free(xsh_resid_tab **resid)
Free memory associated to a resid_tab.
xsh_resid_tab * xsh_resid_tab_load(cpl_frame *resid_tab_frame)
Load a residual tab from a frame.
cpl_vector * xsh_spectralformat_list_get_orders(xsh_spectralformat_list *list, float lambda)
Returns list of absolute orders containing lambda.
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.
void xsh_the_map_free(xsh_the_map **list)
free memory associated to a the_map
void xsh_the_map_set_arcline(xsh_the_map *list, int idx, float wavelength, int order, int slit_index, float slit_position, double detx, double dety)
xsh_the_map * xsh_the_map_create(int size)
Create an empty theoretical map.
cpl_frame * xsh_the_map_save(xsh_the_map *list, const char *filename)
save a the_map to a frame
void xsh_wavesol_apply_shift(xsh_wavesol *wsol, float shift_x, float shift_y)
Apply a shift on X and Y to wave solution.
cpl_table * xsh_wavesol_trace(xsh_wavesol *wsol, double *lambda, double *order, double *slit, int size)
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
cpl_frame * xsh_wavesol_save(xsh_wavesol *w, cpl_table *trace, const char *filename, const char *tag)
save a wavelength solution
#define XSH_REGDEBUG(...)
Definition: xsh_error.h:132
#define check(COMMAND)
Definition: xsh_error.h:71
#define XSH_ASSURE_NOT_NULL(pointer)
Definition: xsh_error.h:99
cpl_frame * xsh_afcthetab_create(cpl_frame *wave_tab_frame, cpl_frame *model_config_frame, int order, cpl_frame *spectralformat_frame, cpl_frame *arclines_frame, int xmin, int ymin, xsh_instrument *instr, const int clean_tmp)
Create a The tab for AFC.
Definition: xsh_flexcor.c:363
cpl_frame * xsh_flexcor(cpl_frame *afc_frame, cpl_frame *wave_tab_frame, cpl_frame *model_config_frame, cpl_frame *order_tab_frame, cpl_frame *attresidtab_frame, int afc_xmin, int afc_ymin, xsh_instrument *instr, cpl_frame **afc_order_tab_frame, cpl_frame **afc_model_config_frame)
This function applies the computed shift betwwen AFC CAL and AFC ATT frame.
Definition: xsh_flexcor.c:91
const char * xsh_instrument_arm_tostring(xsh_instrument *i)
Get the string associated with an arm.
XSH_INSTRCONFIG * xsh_instrument_get_config(xsh_instrument *i)
Get the instrument default set of keywords.
XSH_ARM xsh_instrument_get_arm(xsh_instrument *i)
Get an arm 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
cpl_table * xsh_model_io_output_cfg(struct xs_3 *p_xs_3)
int xsh_model_offset(DOUBLE slit_pix_shift, DOUBLE disp_pix_shift, struct xs_3 *p_xs_3)
convert a pixel shift measured on the detector to a shift in detector centroid
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.
cpl_frame * xsh_model_pipe_anneal(cpl_frame *cfg_frame, cpl_frame *resid_frame, int maxit, double ann_fac, int scenario, int rec_id)
Run the annealing (optimisation) algoritm to improve the fit of the model parameter set to a given wa...
int size
int * y
int * x
#define xsh_msg(...)
Print a message on info level.
Definition: xsh_msg.h:121
void xsh_pfits_set_pcatg(cpl_propertylist *plist, const char *value)
Write the PCATG value.
Definition: xsh_pfits.c:1008
void xsh_free_vector(cpl_vector **v)
Deallocate a vector and set the pointer to NULL.
Definition: xsh_utils.c:2284
void xsh_free_frame(cpl_frame **f)
Deallocate a frame and set the pointer to NULL.
Definition: xsh_utils.c:2269
void xsh_free_table(cpl_table **t)
Deallocate a table and set the pointer to NULL.
Definition: xsh_utils.c:2133
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_ARM_UNDEFINED
@ XSH_ARM_NIR
int xsh_resid_tab_get_size(xsh_resid_tab *resid)
int order
Definition: xsh_detmon_lg.c:80
#define XSH_MOD_CFG_OPT_AFC
Definition: xsh_dfs.h:1256
#define XSH_COL_MODEL_CONF_BEST
Definition: xsh_dfs.h:1294
#define XSH_WAVE_TAB_AFC
Definition: xsh_dfs.h:557
#define XSH_ORDER_TAB_AFC
Definition: xsh_dfs.h:712
#define XSH_GET_TAG_FROM_ARM(TAG, instr)
Definition: xsh_dfs.h:1548
#define XSH_COL_MODEL_CONF_NAME
Definition: xsh_dfs.h:1292
#define XSH_GET_TAG_FROM_MODE(TAG, instr)
Definition: xsh_dfs.h:1594
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_AFC_XSHIFT_C
Definition: xsh_pfits_qc.h:60
#define XSH_QC_AFC_YSHIFT
Definition: xsh_pfits_qc.h:58
#define XSH_QC_AFC_XSHIFT
Definition: xsh_pfits_qc.h:57
#define XSH_QC_AFC_YSHIFT_C
Definition: xsh_pfits_qc.h:61