X-shooter Pipeline Reference Manual 3.8.15
xsh_utils_response.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 * $Author: amodigli $
21 * $Date: 2013-07-12 15:10:08 $
22 * $Revision: 1.34 $
23 */
24
25#ifdef HAVE_CONFIG_H
26# include <config.h>
27#endif
28
29#include <math.h>
30#include <xsh_utils_response.h>
32#include <xsh_dfs.h>
33#include <xsh_pfits_qc.h>
34#include <xsh_error.h>
35/*--------------------------------------------------------------------------*/
41/*--------------------------------------------------------------------------*/
44/*---------------------------------------------------------------------------
45 Includes
46 ---------------------------------------------------------------------------*/
47
48static double
49xsh_resample_double(double wout,double* pw1,double* pf1,double wmin,double wmax,int size_obs)
50{
51
52 double m=0;
53 double w1=0;
54 double w2=0;
55 double f1=0;
56 double f2=0;
57 double fout=0;
58
59 int i=0;
60
61 /* are we in case of extrapolation ? */
62 if( wout < wmin ) {
63 /* w < wmin */
64 m=(pf1[1]-pf1[0])/(pw1[1]-pw1[0]);
65 fout=m*(wout-wmin)+pf1[0];
66 } else if ( wout > wmax ) {
67 /* w > wmax */
68 m=(pf1[size_obs-1]-pf1[size_obs-2])/(pw1[size_obs-1]-pw1[size_obs-2]);
69 fout=m*(wout-wmax)+pf1[size_obs-1];
70
71 } else {
72 /* else, interpolate */
73 for(i=0;i<size_obs;i++) {
74 if(pw1[i]<wout) {
75 w1=pw1[i];
76 f1=pf1[i];
77 } else {
78 w2=pw1[i];
79 f2=pf1[i];
80 break;
81 }
82 }
83 m=(f2-f1)/(w2-w1);
84 fout=m*(wout-w1)+f1;
85 }
86 return fout;
87}
88
89static cpl_table*
90xsh_table_resample_table(cpl_table* tinp,const char* cwinp, const char* cfinp,
91 cpl_table* tref,const char* cwref, const char* cfref)
92{
93
94 cpl_table* tres=NULL;
95 double* pwinp=NULL;
96 double* pwref=NULL;
97 double* pfinp=NULL;
98 double* pfres=NULL;
99 int size_inp=0;
100 int size_ref=0;
101 double wmin=0;
102 double wmax=0;
103 int i=0;
104
105 check(size_inp=cpl_table_get_nrow(tinp));
106 check(size_ref=cpl_table_get_nrow(tref));
107 //xsh_msg("cwinp=%s",cwinp);
108 check(wmin=cpl_table_get_column_min(tinp,cwinp));
109 check(wmax=cpl_table_get_column_max(tinp,cwinp));
110 //xsh_msg("wmin=% vag wmax=%g",wmin,wmax);
111
112 /* we first duplicate ref to resampled table to have the same wavelength
113 sampling and then we take care of computing the flux */
114 tres=cpl_table_duplicate(tref);
115 check(pwinp=cpl_table_get_data_double(tinp,cwinp));
116 check(pwref=cpl_table_get_data_double(tref,cwref));
117
118 check(pfinp=cpl_table_get_data_double(tinp,cfinp));
119 check(pfres=cpl_table_get_data_double(tres,cfref));
120
121 for(i=0;i<size_ref;i++) {
122 pfres[i]=xsh_resample_double(pwref[i],pwinp,pfinp,wmin,wmax,size_inp);
123 }
124
125 cleanup:
126 return tres;
127}
128
129
130
131
132
133static cpl_table*
134xsh_table_downsample_table(cpl_table* tinp,const char* cwinp, const char* cfinp,
135 cpl_table* tref,const char* cwref, const char* cfref)
136{
137
138 cpl_table* tres=NULL;
139 double* pwinp=NULL;
140 double* pwref=NULL;
141 double* pwres=NULL;
142
143 double* pfinp=NULL;
144 //double* pfref=NULL;
145 double* pfres=NULL;
146
147 int size_inp=0;
148 int size_ref=0;
149 double wmin=0;
150 double wmax=0;
151 int i=0;
152
153 check(size_inp=cpl_table_get_nrow(tinp));
154 check(size_ref=cpl_table_get_nrow(tref));
155 //xsh_msg("cwinp=%s",cwinp);
156 check(wmin=cpl_table_get_column_min(tinp,cwinp));
157 check(wmax=cpl_table_get_column_max(tinp,cwinp));
158 //xsh_msg("wmin=%g wmax=%g",wmin,wmax);
159
160 /* we create a table of size as large as the reference one */
161 tres=cpl_table_new(size_ref);
162 /* we create the minimum tables required */
163 cpl_table_new_column(tres,cwref,CPL_TYPE_DOUBLE);
164 cpl_table_new_column(tres,cfref,CPL_TYPE_DOUBLE);
165 cpl_table_fill_column_window_double(tres,cwref,0,size_ref,0.);
166 cpl_table_fill_column_window_double(tres,cfref,0,size_ref,0.);
167
168 //cpl_table_dump_structure(tres,stdout);
169 check(pwinp=cpl_table_get_data_double(tinp,cwinp));
170 check(pwref=cpl_table_get_data_double(tref,cwref));
171
172 check(pfinp=cpl_table_get_data_double(tinp,cfinp));
173 check(pfres=cpl_table_get_data_double(tres,cfref));
174
175 check(pwres=cpl_table_get_data_double(tres,cwref));
176 check(pfres=cpl_table_get_data_double(tres,cfref));
177
178 for(i=0;i<size_ref;i++) {
179 pwres[i]=pwref[i];
180 pfres[i]=xsh_resample_double(pwres[i],pwinp,pfinp,wmin,wmax,size_inp);
181 }
182
183 cleanup:
184 return tres;
185}
186
187//re-sample a table to a uniform sampling step */
188cpl_table*
189xsh_table_resample_uniform(cpl_table* tinp,const char* cwinp, const char* cfinp,
190 const double wstp)
191{
192
193 cpl_table* tres=NULL;
194 double* pwinp=NULL;
195 double* pwref=NULL;
196 double* pfinp=NULL;
197 double* pfres=NULL;
198 int size_inp=0;
199 double wmin=0;
200 double wmax=0;
201 int i=0;
202 int size=0;
203 check(size_inp=cpl_table_get_nrow(tinp));
204 //xsh_msg("cwinp=%s",cwinp);
205 check(wmin=cpl_table_get_column_min(tinp,cwinp));
206 check(wmax=cpl_table_get_column_max(tinp,cwinp));
207
208 size=(int)((wmax-wmin)/wstp+0.5);
209 //xsh_msg("wmin=%15.10g wmax=%15.10g wstp=%15.10g size=%d",wmin,wmax,wstp,size);
210 /* we first duplicate ref to resampled table to have the same wavelength
211 sampling and then we take care of computing the flux */
212 tres=cpl_table_new(size);
213 cpl_table_new_column(tres,cwinp,CPL_TYPE_DOUBLE);
214 cpl_table_new_column(tres,cfinp,CPL_TYPE_DOUBLE);
215
216 cpl_table_fill_column_window_double(tres,cwinp,0,size,0.);
217 cpl_table_fill_column_window_double(tres,cfinp,0,size,0.);
218 check(pwref=cpl_table_get_data_double(tres,cwinp));
219 for(i=0;i<size;i++) {
220 pwref[i]=wmin+i*wstp;
221 }
222 check(pwinp=cpl_table_get_data_double(tinp,cwinp));
223 check(pfinp=cpl_table_get_data_double(tinp,cfinp));
224 check(pfres=cpl_table_get_data_double(tres,cfinp));
225
226 for(i=0;i<size;i++) {
227 pfres[i]=xsh_resample_double(pwref[i],pwinp,pfinp,wmin,wmax,size_inp);
228 }
229
230 cleanup:
231 return tres;
232}
233
234/*
235static cpl_error_code
236xsh_spectrum_free_table(cpl_table* t)
237{
238
239 cpl_table_unwrap(t,"flux");
240 cpl_table_unwrap(t,"errs");
241 cpl_table_unwrap(t,"qual");
242 cpl_table_unwrap(t,"logwave");
243
244 return cpl_error_get_code();
245}
246*/
247
248static cpl_table*
250{
251
252 cpl_table* table_s=NULL;
253 int size=0;
254 int i=0;
255
256 double wmin=0;
257 double wstp=0;
258
259 double* pwave=NULL;
260
261 double* perrs=NULL;
262 double* pflux=NULL;
263 int* pqual=NULL;
264
265 double* errs=NULL;
266 double* flux=NULL;
267 int* qual=NULL;
268
272 //xsh_msg("wmin=%g[nm] wstp=%g[nm]",wmin,wstp);
273
274 //xsh_msg("Lambda size spectrum=%d",size);
275 table_s=cpl_table_new(size);
276 cpl_table_new_column(table_s,"wavelength",CPL_TYPE_DOUBLE);
277 cpl_table_fill_column_window_double(table_s,"wavelength",0,size,0);
278 pwave=cpl_table_get_data_double(table_s,"wavelength");
279
280 cpl_table_new_column(table_s,"flux",CPL_TYPE_DOUBLE);
281 cpl_table_fill_column_window_double(table_s,"flux",0,size,0);
282 pflux=cpl_table_get_data_double(table_s,"flux");
283
284 cpl_table_new_column(table_s,"errs",CPL_TYPE_DOUBLE);
285 cpl_table_fill_column_window_double(table_s,"errs",0,size,0);
286 perrs=cpl_table_get_data_double(table_s,"errs");
287
288 cpl_table_new_column(table_s,"qual",CPL_TYPE_INT);
289 cpl_table_fill_column_window_int(table_s,"qual",0,size,0);
290 pqual=cpl_table_get_data_int(table_s,"qual");
291
292
296
297 for(i=0;i<size;i++) {
298 pwave[i]=wmin+i*wstp;
299 pflux[i]=flux[i];
300 perrs[i]=errs[i];
301 pqual[i]=qual[i];
302 }
303
304 cpl_table_duplicate_column(table_s,"logwave",table_s,"wavelength");
305 cpl_table_logarithm_column(table_s,"logwave",CPL_MATH_E);
306
307 return table_s;
308
309}
310
311static cpl_error_code xsh_evaluate_tell_model(cpl_table* corr,
312 xsh_instrument* instrument, const int ext, double* mean, double* rms) {
313
314 HIGH_ABS_REGION * phigh = NULL;
315 int nsamples = 0;
316 cpl_table* tab_tmp = NULL;
317 cpl_table* tab_stack = NULL;
318 char fname[256];
320 /* we remove from table points that should not be sampled */
321 if (phigh != NULL) {
322 int counter = 0;
323 //int kh = 0;
324 int nrow = 0;
325 //xsh_msg("Extract from correction only sub-regions to compute residuals" );
326 for (; phigh->lambda_min != 0.; phigh++) {
327 //xsh_msg("Flag [%g,%g]",phigh->lambda_min,phigh->lambda_max);
328 nrow = cpl_table_and_selected_double(corr, "wave", CPL_NOT_LESS_THAN,
329 phigh->lambda_min);
330 //xsh_msg("Select min=%d",nrow);
331 nrow = cpl_table_and_selected_double(corr, "wave", CPL_LESS_THAN,
332 phigh->lambda_max);
333 //xsh_msg("Select max=%d",nrow);
334 tab_tmp = cpl_table_extract_selected(corr);
335 cpl_table_select_all(corr);
336 if (counter == 0) {
337 //xsh_msg("create stack table nraws=%d",cpl_table_get_nrow(tab_corr));
338 tab_stack = cpl_table_duplicate(tab_tmp);
339 counter++;
340 } else {
341 nrow = cpl_table_get_nrow(tab_stack);
342 //xsh_msg("append to stack table size=%d",nrow);
343 cpl_table_insert(tab_stack, tab_tmp, nrow);
344 }
345 xsh_free_table(&tab_tmp);
346 nsamples++;
347 }
348
349 }
350 sprintf(fname, "tab_corr_sampl.fits");
351 //xsh_table_save(tab_stack, NULL, NULL, fname, ext);
352
353 *mean = fabs( ( cpl_table_get_column_mean(tab_stack, "correction") - 1. ) );
354 *rms = cpl_table_get_column_stdev(tab_stack, "correction");
355
356 //xsh_msg("mean=%g rms=%g", *mean, *rms);
357 xsh_free_table(&tab_tmp);
358 xsh_free_table(&tab_stack);
359 return cpl_error_get_code();
360
361}
362
363static cpl_table*
364xsh_extract_ranges_to_fit(cpl_table* table_in,const char* cwav,xsh_instrument* instrument) {
365
366 cpl_table* stack = NULL;
367 HIGH_ABS_REGION * phigh = NULL;
368
370
371 cpl_table* tab_tmp = NULL;
372 int nsamples = 0;
373 /*
374 cpl_table_dump(table_in,0,2,stdout);
375 xsh_msg("wmin=%g wmax=%g",
376 cpl_table_get_column_min(table_in,cwav),
377 cpl_table_get_column_max(table_in,cwav));
378 */
379 /* we remove from table points that should not be sampled */
380 if (phigh != NULL) {
381 int counter = 0;
382 //int kh = 0;
383 int nrow = 0;
384 //xsh_msg("Extract from ratio only regions to be sampled" );
385 for (; phigh->lambda_min != 0.; phigh++) {
386 //xsh_msg("Flag [%g,%g]",phigh->lambda_min,phigh->lambda_max);
387 check(nrow = cpl_table_and_selected_double(table_in, cwav, CPL_NOT_LESS_THAN,
388 phigh->lambda_min));
389 //xsh_msg("Select min=%d",nrow);
390 check(nrow = cpl_table_and_selected_double(table_in,cwav, CPL_LESS_THAN,
391 phigh->lambda_max));
392 //xsh_msg("Select max=%d",nrow);
393 check(tab_tmp = cpl_table_extract_selected(table_in));
394 cpl_table_select_all(table_in);
395 if (counter == 0) {
396 //xsh_msg("create stack table");
397 stack = cpl_table_duplicate(tab_tmp);
398 counter++;
399 } else {
400 nrow=cpl_table_get_nrow(stack);
401 //xsh_msg("append to stack table size=%d",nrow);
402 cpl_table_insert(stack, tab_tmp, nrow);
403 }
404 xsh_free_table(&tab_tmp);
405 nsamples++;
406 }
407 }
408 xsh_free_table(&tab_tmp);
409
410 cleanup:
411 return stack;
412}
413/*
414static cpl_error_code
415xsh_extract_points_to_fit(cpl_table* table,const char* cwav,const char* cratio, int nsamples, xsh_instrument* instrument,cpl_vector** vec_wave,cpl_vector** vec_flux)
416{
417
418
419 double* pvec_wave = NULL;
420 double* pvec_flux = NULL;
421 int kh=0;
422 HIGH_ABS_REGION * phigh = NULL;
423 cpl_table* tab_tmp=NULL;
424 double wmin=0;
425 double wmax=0;
426 double wstp=1.;
427 int nsel=0;
428 double flux_sto=0;
429 double wmin_fit_region=9999;
430 double wmax_fit_region=-9999;
431 int kh_min=0;
432 int kh_max=0;
433
434 wmin=cpl_table_get_column_min(table,cwav);
435 wmax=cpl_table_get_column_max(table,cwav);
436
437 phigh=xsh_fill_tell_fit_regions(instrument,NULL);
438 for (kh = 0; phigh->lambda_min != 0.; phigh++) {
439 if( wmin_fit_region> phigh->lambda_min) wmin_fit_region=phigh->lambda_min;
440 if( wmax_fit_region< phigh->lambda_max) wmax_fit_region=phigh->lambda_max;
441 }
442 xsh_msg("wmin_fit_region=%g wmax_fit_region=%g",wmin_fit_region,wmax_fit_region);
443 xsh_msg("wmin=%g wmax=%g",wmin,wmax);
444 //if(wmin_fit_region<wmin) {
445 nsamples+=1;
446 kh_min=1;
447 //}
448
449 if(wmax_fit_region>wmax) {
450 nsamples+=1;
451 kh_max=nsamples-1;
452 }
453
454
455 *vec_wave = cpl_vector_new(nsamples);
456 *vec_flux = cpl_vector_new(nsamples);
457 pvec_wave = cpl_vector_get_data(*vec_wave);
458 pvec_flux = cpl_vector_get_data(*vec_flux);
459 //xsh_msg("sk1 nsamples=%d",nsamples);
460
461 //if(wmin_fit_region<wmin) {
462 // add one point at the start of the wave range
463 cpl_table_and_selected_double(table, cwav, CPL_LESS_THAN,wmin+wstp);
464 tab_tmp = cpl_table_extract_selected(table);
465
466 pvec_wave[0] = wmin;
467 check(pvec_flux[0] = cpl_table_get_column_median(tab_tmp, cratio));
468 cpl_table_select_all(table);
469 xsh_free_table(&tab_tmp);
470 //}
471
472 phigh=xsh_fill_tell_fit_regions(instrument,NULL);
473 if (phigh != NULL) {
474 //int nrow = 0;
475 //xsh_msg("Fill median ratio values" );
476 for (kh = kh_min; phigh->lambda_min != 0.; phigh++) {
477 //xsh_msg("Process line %d Flag [%g,%g]",kh,phigh->lambda_min,phigh->lambda_max);
478 cpl_table_and_selected_double(table, cwav,
479 CPL_NOT_LESS_THAN, phigh->lambda_min);
480 //xsh_msg("Select min=%d",nrow);
481 nsel=cpl_table_and_selected_double(table, cwav, CPL_LESS_THAN,
482 phigh->lambda_max);
483 //xsh_msg("Select max=%d",nrow);
484 if(nsel>0) {
485 tab_tmp = cpl_table_extract_selected(table);
486 cpl_table_select_all(table);
487
488 pvec_wave[kh] = 0.5 * (phigh->lambda_max + phigh->lambda_min);
489 check(pvec_flux[kh] = cpl_table_get_column_median(tab_tmp, cratio));
490 flux_sto=pvec_flux[kh];
491 xsh_free_table(&tab_tmp);
492 //xsh_msg("wave=%g",pvec_wave[kh+1]);
493 } else {
494 pvec_wave[kh] = 0.5 * (phigh->lambda_max + phigh->lambda_min);
495 pvec_flux[kh]=flux_sto;
496 //xsh_msg("wave=%g",pvec_wave[kh+1]);
497 }
498 kh++;
499 //TO_BE_FIXED
500 }
501 }
502
503 if(wmax_fit_region>wmax) {
504 cpl_table_select_all(table);
505
506 cpl_table_and_selected_double(table, cwav, CPL_NOT_LESS_THAN,wmax-wstp);
507 tab_tmp = cpl_table_extract_selected(table);
508
509 pvec_wave[kh_max] = wmax;
510 xsh_msg("wave=%g",pvec_wave[kh+1]);
511 check(pvec_flux[kh_max] = cpl_table_get_column_median(tab_tmp, cratio));
512 xsh_free_table(&tab_tmp);
513 }
514
515 cleanup:
516
517 return cpl_error_get_code();
518}
519*/
520
521static cpl_error_code
522xsh_extract_points_to_fit(cpl_table* table,const char* cwav,const char* cratio, const int nsamples, xsh_instrument* instrument,cpl_vector** vec_wave,cpl_vector** vec_flux)
523{
524
525
526 double* pvec_wave = NULL;
527 double* pvec_flux = NULL;
528 int kh=0;
529 HIGH_ABS_REGION * phigh = NULL;
530 cpl_table* tab_tmp=NULL;
531 double wmin=0;
532 double wmax=0;
533 double wstp=1.;
534 /*
535 double wmin_fit_region=9999;
536 double wmax_fit_region=-9999;
537 */
538
539 wmin=cpl_table_get_column_min(table,cwav);
540 wmax=cpl_table_get_column_max(table,cwav);
541
542 *vec_wave = cpl_vector_new(nsamples+2);
543 *vec_flux = cpl_vector_new(nsamples+2);
544
545 /*
546 phigh=xsh_fill_tell_fit_regions(instrument,NULL);
547 for (kh = 0; phigh->lambda_min != 0.; phigh++) {
548 if( wmin_fit_region> phigh->lambda_min) wmin_fit_region=phigh->lambda_min;
549 if( wmax_fit_region< phigh->lambda_max) wmax_fit_region=phigh->lambda_max;
550 }
551
552 if(wmax_fit_region>wmax) {
553 *vec_wave = cpl_vector_new(nsamples+2);
554 *vec_flux = cpl_vector_new(nsamples+2);
555 } else {
556 *vec_wave = cpl_vector_new(nsamples+1);
557 *vec_flux = cpl_vector_new(nsamples+1);
558 }
559 */
560 pvec_wave = cpl_vector_get_data(*vec_wave);
561 pvec_flux = cpl_vector_get_data(*vec_flux);
562 //xsh_msg("sk1 nsamples=%d",nsamples);
563
564 // add one point at the start of the wave range
565 cpl_table_and_selected_double(table, cwav, CPL_LESS_THAN,wmin+wstp);
566 tab_tmp = cpl_table_extract_selected(table);
567 pvec_wave[0] = wmin;
568 pvec_flux[0] = cpl_table_get_column_median(tab_tmp, cratio);
569 cpl_table_select_all(table);
570 xsh_free_table(&tab_tmp);
571
573 if (phigh != NULL) {
574 //int nrow = 0;
575 //xsh_msg("Fill median ratio values" );
576 for (kh = 0; phigh->lambda_min != 0.; phigh++) {
577 //xsh_msg("Process line %d Flag [%g,%g]",kh,phigh->lambda_min,phigh->lambda_max);
578 cpl_table_and_selected_double(table, cwav,
579 CPL_NOT_LESS_THAN, phigh->lambda_min);
580 //xsh_msg("Select min=%d",nrow);
581 cpl_table_and_selected_double(table, cwav, CPL_LESS_THAN,
582 phigh->lambda_max);
583 //xsh_msg("Select max=%d",nrow);
584 tab_tmp = cpl_table_extract_selected(table);
585 cpl_table_select_all(table);
586
587 pvec_wave[kh+1] = 0.5 * (phigh->lambda_max + phigh->lambda_min);
588 pvec_flux[kh+1] = cpl_table_get_column_median(tab_tmp, cratio);
589 kh++;
590 xsh_free_table(&tab_tmp);
591 }
592 }
593
594 /*
595 if(wmax_fit_region>wmax) {
596 cpl_table_select_all(table);
597
598 cpl_table_and_selected_double(table, cwav, CPL_NOT_LESS_THAN,wmax-wstp);
599 tab_tmp = cpl_table_extract_selected(table);
600
601 pvec_wave[kh+1] = wmax;
602 xsh_msg("wave=%g",pvec_wave[kh+1]);
603 check(pvec_flux[kh+1] = cpl_table_get_column_median(tab_tmp, cratio));
604 xsh_free_table(&tab_tmp);
605 }
606 */
607
608
609 cpl_table_select_all(table);
610
611 cpl_table_and_selected_double(table, cwav, CPL_NOT_LESS_THAN,wmax-wstp);
612 tab_tmp = cpl_table_extract_selected(table);
613
614 pvec_wave[kh+1] = wmax;
615 pvec_flux[kh+1] = cpl_table_get_column_median(tab_tmp, cratio);
616 xsh_free_table(&tab_tmp);
617
618
619 //cleanup:
620
621 return cpl_error_get_code();
622}
623
624static cpl_error_code
625xsh_get_xcorrel_peak(cpl_vector* wcorr, cpl_vector* fcorr,XSH_GAUSSIAN_FIT* gfit,const double range,const int ext)
626{
627 int size_x;
628 cpl_table* tab = NULL;
629 cpl_table* ext_tab;
630 char fname[256];
631 cpl_vector* wave=NULL;
632 cpl_vector* flux=NULL;
633
634 size_x=cpl_vector_get_size(wcorr);
635 tab = cpl_table_new(size_x);
636 cpl_table_wrap_double(tab, cpl_vector_get_data(wcorr), "logwave");
637 cpl_table_wrap_double(tab, cpl_vector_get_data(fcorr), "flux");
638
639 sprintf(fname, "fcorr_org.fits");
640 //xsh_table_save(tab,NULL,NULL,fname,ext);
641 cpl_table_and_selected_double(tab, "logwave", CPL_GREATER_THAN,
642 (*gfit).peakpos - range);
643 cpl_table_and_selected_double(tab, "logwave", CPL_LESS_THAN,
644 (*gfit).peakpos + range);
645
646 ext_tab = cpl_table_extract_selected(tab);
647 cpl_table_unwrap(tab,"logwave");
648 cpl_table_unwrap(tab,"flux");
649 xsh_free_table(&tab);
650 sprintf(fname, "fcorr_ext.fits");
651 //xsh_table_save(ext_tab,NULL,NULL,fname,ext);
652 int next = 0;
653 next = cpl_table_get_nrow(ext_tab);
654 sprintf(fname, "fcorr_tab.fits");
655 //xsh_table_save(ext_tab,NULL,NULL,fname,ext);
656 double peakpos=0;
657 double sigma=0;
658 double area=0;
659 double offset=0;
660 double mse=0;
661
662 wave = cpl_vector_wrap(next, cpl_table_get_data_double(ext_tab, "logwave"));
663 flux = cpl_vector_wrap(next, cpl_table_get_data_double(ext_tab, "flux"));
664 //xsh_msg("gauss before fit: peak=%g sigma=%g",(*gfit).peakpos,(*gfit).sigma);
665
666 cpl_vector_fit_gaussian( wave, NULL, flux,NULL, CPL_FIT_ALL,
667 &peakpos, &sigma,&area, &offset, &mse,NULL,NULL);
668 //xsh_msg("gauss check fit: peak=%g sigma=%g",peakpos,sigma);
669 cpl_vector_fit_gaussian( wave, NULL, flux,NULL, CPL_FIT_ALL,
670 &gfit->peakpos, &gfit->sigma,&gfit->area, &gfit->offset, &gfit->mse,NULL,NULL);
671 //double fwhm = CPL_MATH_FWHM_SIG * gfit.sigma;
672 //xsh_msg("gauss after fit: peak=%g sigma=%g",(*gfit).peakpos,(*gfit).sigma);
673
674 //cleanup:
675 cpl_vector_unwrap(wave);
676 cpl_vector_unwrap(flux);
677 xsh_free_table(&ext_tab);
678 return cpl_error_get_code();
679}
680
681static cpl_table*
682xsh_table_select_range(cpl_table* table_in,const char* col,const double wmin, const double wmax)
683{
684 cpl_table* table_ext=NULL;
685
686 cpl_table_and_selected_double(table_in, col, CPL_NOT_LESS_THAN, wmin);
687 cpl_table_and_selected_double(table_in, col, CPL_NOT_GREATER_THAN, wmax);
688 table_ext=cpl_table_extract_selected(table_in);
689
690 return table_ext;
691}
692
693/* perform correlation of spectrum and model over a given range. Correlation peak determination is refined by Gauss fit
694 * result is the shift to be applied and the FWHM of the correlation function (useful to make 2nd iteration of correlation
695 * near peak pos)
696 */
697cpl_error_code
698xsh_correl_spectra(double* flux_s, double* flux_m, const int size,
699 const int hsearch, const double wlogstp,const int norm_xcorr,
700 const double range, const int ext,XSH_GAUSSIAN_FIT* gfit)
701{
702 double* xcorr = NULL;
703 double* pwcorr = NULL;
704 double corr = 0;
705 double delta=0;
706 char fname[256];
707 int i=0;
708 double size_corr = 2 * hsearch + 1;
709 cpl_vector* f1=NULL;
710 cpl_vector* f2=NULL;
711 cpl_vector* correl=NULL;
712 double shift=0;
713 cpl_vector* fcorr = NULL;
714 cpl_vector* wcorr = NULL;
715
716 check(xcorr = xsh_function1d_xcorrelate(flux_s, size, flux_m, size, hsearch,
717 norm_xcorr,&corr, &delta));
718
719
720 check(f1=cpl_vector_wrap(size, flux_s));
721 f2=cpl_vector_wrap(size, flux_m);
722 //cpl_vector_save(f1,"f1.fits",CPL_BPP_IEEE_FLOAT,NULL,CPL_IO_DEFAULT);
723 //cpl_vector_save(f2,"f2.fits",CPL_BPP_IEEE_FLOAT,NULL,CPL_IO_DEFAULT);
724 //xsh_msg("size=%d",size);
725 correl=cpl_vector_new(size_corr);
726 check(shift=cpl_vector_correlate(correl,f1,f2));
727 cpl_vector_unwrap(f1);
728 cpl_vector_unwrap(f2);
729 //cpl_vector_save(correl,"correl.fits",CPL_BPP_IEEE_FLOAT,NULL,CPL_IO_DEFAULT);
730 xsh_msg("shift=%g",shift);
731 //xsh_msg("correlation value=%g measured shift[sample units]=%g measured shift[lognm]=%g", corr, delta,delta*);
732 //xsh_msg("hsearch=%d delta=%g sampling step=%g",hsearch,delta,wlogstp);
733
734 gfit->peakpos = (hsearch + delta)*wlogstp;
735 gfit->sigma = wlogstp*10;
736 gfit->area = 1.;
737
738 xsh_msg("gauss guess: peak: %12.8g sigma %g", gfit->peakpos, gfit->sigma);
739
740 check(fcorr = cpl_vector_wrap(size_corr, xcorr));
741 /* AMO changed: to remove a leak (next we override values and
742 * later we do cpl_vector_unwrap
743 */
744 pwcorr= (double*) cpl_calloc(size_corr, sizeof(double));
745 for (i = 0; i < size_corr; i++) {
746 pwcorr[i] = i*wlogstp;
747 }
748 check(wcorr = cpl_vector_wrap(size_corr, pwcorr));
749 sprintf(fname,"fcorr.fits");
750 //xsh_vector_save(fcorr,NULL,fname,ext);
751 sprintf(fname,"wcorr.fits");
752 //check(xsh_vector_save(wcorr,NULL,fname,ext));
753
754 /* 4) fit Gaussian to +/-0.001 in log_lam_nm around cross correlation peak */
755 check(xsh_get_xcorrel_peak(wcorr,fcorr,gfit,range,ext));
756
757 xsh_msg("gauss fit: peak[lognm]: %12.8g sigma[lognm] %g peak[sampl_units]: %12.8g sigma[sampl_units] %g",
758 gfit->peakpos, gfit->sigma,gfit->peakpos/wlogstp, gfit->sigma/wlogstp);
759
760
761 cleanup:
762 cpl_vector_unwrap(fcorr);
763 cpl_vector_unwrap(wcorr);
764 xsh_free_vector(&correl);
765 cpl_free(xcorr);
766 cpl_free(pwcorr);
767 return cpl_error_get_code();
768
769}
770static cpl_error_code
771xsh_align_model_to_spectrum(cpl_table* table_me,cpl_table* table_se,const double wlogstp,const int norm_xcorr,const int ext,cpl_table** table_mm)
772{
773 cpl_table* table_rm = NULL;
774 cpl_table* table_rs = NULL;
775 char fname[256];
776 int i=0;
777 int size_mm=cpl_table_get_nrow(*table_mm);
778 //cpl_lowpass filter_type = CPL_LOWPASS_LINEAR;
779 cpl_vector* vfluxm = NULL;
780 cpl_vector* vlowpass_fluxm = NULL;
781
782 /* TODO: this re-sampling was decided not to use */
783 //check(cpl_table_dump_structure(table_me,stdout));
784 //xsh_msg("logspt=%g",wlogstp);
785 check(table_rm = xsh_table_resample_uniform(table_me, "logwave", "flux", wlogstp));
786 //xsh_msg("size of spectrum tab=%d",(int)cpl_table_get_nrow(table_rm));
787 sprintf(fname, "model_selected_wrange_log_resampled_uniform.fits");
788 //check(xsh_table_save(table_rm,NULL,NULL,fname,ext));
789
790
791 check(table_rs = xsh_table_downsample_table(table_se, "logwave", "flux", table_rm,
792 "logwave", "flux"));
793
794 sprintf(fname, "spectrum_resampled_used_in_correlation.fits");
795 //xsh_table_save(table_rs,NULL,NULL,fname,ext);
796
797 /* perform correlation of spectrum and model over a given range. Correlation peak determination is refined by Gauss fit
798 * result is the shift to be applied and the FWHM of the correlation function (useful to make 2nd iteration of correlation
799 * near peak pos)
800 */
801
802 //double corr = 0;
803 double* flux_s = NULL;
804 double* flux_m = NULL;
805 int size_s = 0;
806 int hsearch = 200;
807 /*
808 int size_m = 0;
809
810 double delta = 0;
811 double size_x = 2 * hsearch + 1;
812 double* xcorr = NULL;
813 double* pwcorr = NULL;
814 int size=0;
815 */
816 XSH_GAUSSIAN_FIT gfit;
817
818 double wrange=0.0005;
819 //size=cpl_table_get_nrow(table_rs);
820
821 flux_s = cpl_table_get_data_double(table_rs, "flux");
822 flux_m = cpl_table_get_data_double(table_rm, "flux");
823 size_s = cpl_table_get_nrow(table_rs);
824 //size_m = cpl_table_get_nrow(table_me);
825
826 xsh_msg("hsearch=%d",hsearch);
827 check(xsh_correl_spectra(flux_s,flux_m,size_s,hsearch,wlogstp,norm_xcorr,wrange,ext,&gfit));
828 //xsh_msg("sigma=%g",gfit.sigma);
829 xsh_msg("computed shift[lognm]=%g",(gfit.peakpos-hsearch*wlogstp));
830 hsearch=(int)3.*CPL_MATH_FWHM_SIG * gfit.sigma/wlogstp;
831 xsh_msg("hsearch=%d",hsearch);
832 check(xsh_correl_spectra(flux_s,flux_m,size_s,hsearch,wlogstp,norm_xcorr,wrange,ext,&gfit));
833
834 /* 5: adjust model wavelength scale:
835 :log_lam_nm_0 = :log_lam_nm+<position>
836 */
837 cpl_table_add_scalar(*table_mm, "logwave", (gfit.peakpos-hsearch*wlogstp));
838 xsh_msg("computed shift[lognm]=%g",(gfit.peakpos-hsearch*wlogstp));
839 //xsh_table_save(*table_mm,NULL,NULL,"model_shift.fits",ext);
840
841
842 /*
843 6) create Gaussian with measured FWHM (corrected by FWHM of model spectrum,
844 but that effect is small since the model spectra have a resolution of 60000,
845 so the FWHM of the cross-correlation peak is dominated by the observation)
846
847 5) convolve model with Gaussain in log_lam_nm space (X-shooter has contant lambda/dlambda)
848
849 We apply this operations to the model sub-range corresponding to the spectrum
850 */
851
852 /* NB:
853 * the only eventually allowed filer in cpl_vector_filter_lowpass_create is CPL_LOWPASS_LINEAR
854 * as CPL_LOWPASS_GAUSSIAN uses a fixed FWHM. Checking result we later preferred to use
855 * cpl_wlcalib_xc_convolve_create_kernel that is deprecated but gives better results.
856 * There the FWHM parameter is actually the sigma: FWHM/CPL_MATH_FWHM_SIG
857 */
858
859 double* pvlowpass_fluxm = NULL;
860 double* pflux_c = NULL;
861 double fwhm = CPL_MATH_FWHM_SIG * gfit.sigma;
862 double fwhm_nm = pow(CPL_MATH_E, fwhm);
863 cpl_vector * conv_kernel =NULL;
864 //cpl_vector * spec_clean=NULL;
865 //int fwhm_pix=(int)(fwhm_nm+0.50000);
866 int fwhm_pix = (int) (fwhm / wlogstp + 0.5);
867
868 cpl_table_new_column(*table_mm, "flux_c", CPL_TYPE_DOUBLE);
869 cpl_table_fill_column_window_double(*table_mm, "flux_c", 0, size_mm, 0);
870 xsh_msg("fwhm_nm=%12.10g fwhm_pix=%d",fwhm_nm,fwhm_pix);
871
872 vfluxm = cpl_vector_wrap(size_mm,cpl_table_get_data_double(*table_mm, "flux"));
873
874 /* we use CPL deprecated function but it works properly */
875 conv_kernel = cpl_wlcalib_xc_convolve_create_kernel(fwhm_pix/CPL_MATH_FWHM_SIG,fwhm_pix/CPL_MATH_FWHM_SIG);
876 vlowpass_fluxm=cpl_vector_duplicate(vfluxm);
877 cpl_wlcalib_xc_convolve(vlowpass_fluxm, conv_kernel);
878
879 //vlowpass_fluxm = cpl_vector_filter_lowpass_create(vfluxm, filter_type,fwhm_pix);
880
881 pvlowpass_fluxm = cpl_vector_get_data(vlowpass_fluxm);
882 cpl_vector_unwrap(vfluxm);
883
884 pflux_c = cpl_table_get_data_double(*table_mm, "flux_c");
885 for (i = 0; i < size_mm; i++) {
886 pflux_c[i] = pvlowpass_fluxm[i];
887 }
888 sprintf(fname, "model_convolved.fits");
889 //xsh_table_save(*table_mm,NULL,NULL,fname,ext);
890
891 cleanup:
892
893 xsh_free_vector(&conv_kernel);
894 xsh_free_vector(&vlowpass_fluxm);
895 xsh_free_table(&table_rs);
896 xsh_free_table(&table_rm);
897
898 return cpl_error_get_code();
899
900}
901
902cpl_table*
903xsh_telluric_model_eval(cpl_frame* frame_m,xsh_spectrum* s,xsh_instrument* instrument,cpl_size* model_idx)
904{
905
906 cpl_table* tab_res=NULL;
907 cpl_table* table_s=NULL;
908 cpl_table* table_m=NULL;
909 cpl_table* table_se=NULL;
910 cpl_table* table_mm=NULL;
911 cpl_table* table_rm = NULL;
912 cpl_table* table_rs_div_mc=NULL;
913 cpl_table* table_me=NULL;
914 cpl_table* table_rs=NULL;
915 const char* name_model=NULL;
916
917 double wmin=0;
918 double wmax=0;
919 //int i=0;
920 double lmin=7.0; //1096.633 nm
921 double lmax=7.1; //1211.967 nm
922 double wlogstp = 10.e-6;
923
924 double um2nm=1000.; // 1nm
925
926 int ext=0;
927 int next=0;
928 double mean=0;
929 double rms=0;
930 //int kh=0;
931 int nsamples=0;
932 char fname[256];
933 cpl_table* tab_corr = NULL;
934 cpl_table* tab_stack = NULL;
935 //int size_mm=0;
936 cpl_table* ext_tab = NULL;
937 double* pmean=NULL;
938 double* prms=NULL;
939 int* pext=NULL;
940 int norm_xcorr=0;
941 HIGH_ABS_REGION * phigh = NULL;
942 cpl_propertylist* qc_head=NULL;
943
945 lmin=6.828; //(923.342 nm)
946 lmax=6.894; //(986.339 nm)
947 //lmax=6.888; //(986.339 nm)
948 wlogstp = 7.826e-6; //(1 nm)
949 norm_xcorr=1;
950 }
951
954 //xsh_msg("wmin=%g wmax=%g",wmin,wmax);
955 //wstp= xsh_spectrum_get_lambda_step(s);
956 check(table_s=xsh_spectrum_to_table(s));
957 sprintf(fname,"spectrum_obs.fits");
958 //check(xsh_table_save(table_s,NULL,NULL,fname,0));
960
961 for (; phigh->lambda_min != 0.; phigh++) {
962 nsamples++;
963 }
964
965 check(name_model=cpl_frame_get_filename(frame_m));
966 //cpl_frame_set_type(frame_m,CPL_FRAME_TYPE_TABLE);
967 check(next=cpl_frame_get_nextensions(frame_m));
968 check(tab_res=cpl_table_new(next));
969 cpl_table_new_column(tab_res,"id",CPL_TYPE_INT);
970 cpl_table_new_column(tab_res,"mean",CPL_TYPE_DOUBLE);
971 cpl_table_new_column(tab_res,"rms",CPL_TYPE_DOUBLE);
972
973 cpl_table_fill_column_window_int(tab_res, "id", 0, next, 0);
974 cpl_table_fill_column_window_double(tab_res, "mean", 0, next, 0);
975 cpl_table_fill_column_window_double(tab_res, "rms", 0, next, 0);
976 pext=cpl_table_get_data_int(tab_res,"id");
977 pmean=cpl_table_get_data_double(tab_res,"mean");
978 prms=cpl_table_get_data_double(tab_res,"rms");
979 check(next=cpl_table_get_nrow(tab_res));
980
981 //xsh_msg("name_model %s",name_model);
982 for (ext = 0; ext < next; ext++) {
983
984 /* load input model spectrum from appropriate extension */
985 check(table_m = cpl_table_load(name_model, ext+1, 1));
986 sprintf(fname,"model_arm.fits");
987 //check(xsh_table_save(table_m,NULL,NULL,fname,ext));
988 /* the model wavelength are in um units, we need to convert to nm */
989 cpl_table_multiply_scalar(table_m, "lam", um2nm);
990 cpl_table_duplicate_column(table_m, "logwave", table_m, "lam");
991 cpl_table_logarithm_column(table_m, "logwave", CPL_MATH_E);
992 //check(xsh_table_save(table_m,NULL,NULL,fname,ext));
993 /* we restrict model to range covered by observation */
994 table_mm=xsh_table_select_range(table_m,"lam",wmin, wmax);
995 xsh_free_table(&table_m);
996 //size_mm = cpl_table_get_nrow(table_mm);
997 sprintf(fname,"model_range_arm.fits");
998 //check(xsh_table_save(table_mm,NULL,NULL,fname,ext));
999
1000 /* Step 2: extract the proper ranges for cross correlation */
1001 /* Apply this to the model */
1002 table_me=xsh_table_select_range(table_mm,"logwave",lmin, lmax);
1003
1004
1005 /* make a test to verify that exponential of log comes back the same
1006 as original wave column value
1007 check(cpl_table_duplicate_column(table_me,"wav",table_me,"logwave"));
1008 check(cpl_table_exponential_column(table_me,"wav",CPL_MATH_E));
1009 */
1010 sprintf(fname,"model_selected_wrange_log.fits");
1011 //check(xsh_table_save(table_me,NULL,NULL,fname,ext));
1012
1013 /* Apply this to the spectrum */
1014 table_se=xsh_table_select_range(table_s,"logwave",lmin, lmax);
1015 cpl_table_duplicate_column(table_se, "wav", table_se, "logwave");
1016 cpl_table_exponential_column(table_se, "wav", CPL_MATH_E);
1017 sprintf(fname,"spectrum_selected_wrange_log.fits");
1018 //check(xsh_table_save(table_se,NULL,NULL,fname,ext));
1019
1020 /* Step 3: now we need to do the cross-correlation:
1021 Before doing this we need to be sure that each spectrum is sampled
1022 to the same point. We use the model as reference and re-sample (linearly)
1023 the observed spectrum so that we can then compute the cross-correlation
1024
1025 Step 4: fit Gaussian to +/-0.001 in log_lam_nm around cross correlation peak
1026
1027 Step 5: adjust model wavelength scale:
1028 :log_lam_nm_0 = :log_lam_nm+<position>
1029
1030
1031 6) create Gaussian with measured FWHM (corrected by FWHM of model spectrum,
1032 but that effect is small since the model spectra have a resolution of 60000,
1033 so the FWHM of the cross-correlation peak is dominated by the observation)
1034
1035 5) convolve model with Gaussain in log_lam_nm space (X-shooter has contant lambda/dlambda)
1036
1037 We apply this operations to the model sub-range corresponding to the spectrum
1038 */
1039 //xsh_msg("size of model tab=%d",(int)cpl_table_get_nrow(table_me));
1040 //xsh_msg("size of spectrum tab=%d",(int)cpl_table_get_nrow(table_se));
1041 //check(xsh_table_save(table_me,NULL,NULL,"table_me.fits",ext));
1042 //check(xsh_table_save(table_se,NULL,NULL,"table_se.fits",ext));
1043 check(xsh_align_model_to_spectrum(table_me,table_se,wlogstp,norm_xcorr,ext,&table_mm));
1044 xsh_free_table(&table_me);
1045 xsh_free_table(&table_se);
1046
1047 //check(xsh_table_save(table_mm,NULL,NULL,"model_shift.fits",ext));
1048
1049 /*
1050 6) convert convolved and sifted model to lam_nm space
1051 */
1052 check(cpl_table_duplicate_column(table_mm,"wave",table_mm,"logwave"));
1053 check(cpl_table_exponential_column(table_mm,"wave",CPL_MATH_E));
1054 sprintf(fname,"model_aligned_to_obs.fits");
1055
1056 //check(xsh_table_save(table_mm,NULL,NULL,fname,ext));
1057
1058 /*
1059 7) divide observed spectrum by convolved and shifted model
1060 Here we get a telluric-corrected observed spectrum
1061 */
1062 xsh_free_table(&table_rs);
1063
1064 xsh_free_table(&table_rm);
1065 /* In order to divide spectrum by model we need to put them on same sampling steps:
1066 * If we re-sample the model to the coarser spectrum step we introduce artifacts in the re-sampled model spectrum
1067 * Better is to up-sample the observed spectrum to the model. This simply artificially refine the observed spectrum
1068 */
1069
1070 /* NOTE THAT WE NEED table_rpfitm at the very end of this function when we downsample table_rs */
1071 check(table_rm = xsh_table_resample_table(table_mm,"wave","flux",table_s,"wavelength","flux"));
1072 sprintf(fname,"table_rm.fits");
1073 //check(xsh_table_save(table_rm,NULL,NULL,fname,ext));
1074
1075
1076 check(
1077 table_rs=xsh_table_resample_table(table_s,"wavelength","flux",table_mm,"wave","flux"));
1078
1079 //PIPPO
1080 sprintf(fname,"spectrum_resampled_to_model.fits");
1081 //check(xsh_table_save(table_rs,NULL,NULL,fname,ext));
1082
1083 cpl_table_duplicate_column(table_rs, "flux_model_convolved", table_mm, "flux_c");
1084 cpl_table_duplicate_column(table_rs, "ratio", table_rs, "flux");
1085 xsh_free_table(&table_mm);
1086 check(cpl_table_divide_columns(table_rs,"ratio","flux_model_convolved"));
1087 /* NB: the ratio spectrum corresponding to the model that minimizes
1088 the resisiduals is actually the final product we need to use for
1089 the final response computation */
1090
1091 sprintf(fname,"spectrum_observed_divided_by_model_convolved.fits");
1092 check(xsh_table_save(table_rs,NULL,NULL,fname,ext));
1093 table_rs_div_mc=cpl_table_duplicate(table_rs);
1094 /*
1095 8) use predefined continuum points (avoiding regions of strong
1096 telluric absorption as well as known stellar lines) and fit spline
1097 without smoothing
1098
1099 (points.list, email from Sabine).
1100 Take around each of them +/- 1nm.
1101 take the median in the range [-1nm+wc,wc+1nm]
1102 and then fit a spline (order 3: cubic)
1103 No smoothing
1104 This gives the continuum fit
1105 */
1106
1107 check(tab_stack=xsh_extract_ranges_to_fit(table_rs,"wave",instrument));
1108 sprintf(fname,"tab_stack.fits");
1109
1110 //check(xsh_table_save(tab_stack,NULL,NULL,fname,ext));
1111
1112 /* here we make the spline fit of median flux computed in each range
1113 * (to rule out that oscillations in each range creates problems)
1114 * */
1115
1116 cpl_vector* vec_wave = NULL;
1117 cpl_vector* vec_flux = NULL;
1118 //double* pwav=NULL;
1119
1120 //check(xsh_table_save(table_rs,NULL,NULL,"pippo.fits",ext));
1121 check(xsh_extract_points_to_fit(table_rs,"wave","ratio",nsamples,instrument,&vec_wave,&vec_flux));
1122
1123 //pwav = cpl_table_get_data_double(table_rm, "wave");
1124 double* sfit=NULL;
1125 double* pfit=NULL;
1126 //double* pwave=NULL;
1127 //cpl_vector_dump(vec_wave,stdout);
1128 //cpl_vector_dump(vec_flux,stdout);
1129
1130 /*
1131
1132 double* pflux=NULL;
1133
1134 sprintf(fname,"ratio_sel.fits");
1135 check(xsh_table_save(tab_stack,NULL,NULL,fname,ext));
1136 */
1137 int nsel = 0;
1138
1139 check(nsel=cpl_table_get_nrow(tab_stack));
1140 cpl_table_new_column(tab_stack, "fit", CPL_TYPE_DOUBLE);
1141 cpl_table_fill_column_window_double(tab_stack, "fit", 0, nsel, 0);
1142 pfit = cpl_table_get_data_double(tab_stack, "fit");
1143 //pwave = cpl_table_get_data_double(tab_stack, "wave");
1144
1145 /*
1146
1147
1148 pwav = cpl_table_get_data_double(tab_stack, "wave");
1149 pflux = cpl_table_get_data_double(tab_stack, "ratio");
1150 vflux=cpl_vector_wrap(size_mm,pflux);
1151 vflux_med=cpl_vector_filter_median_create(vflux,15);
1152 xsh_msg("size_mm=%d",size_mm);
1153 sprintf(fname,"flux_med.fits");
1154 cpl_vector_save(vflux_med,fname,CPL_BPP_IEEE_FLOAT,NULL,CPL_IO_DEFAULT);
1155 mflux=cpl_vector_get_data(vflux_med);
1156
1157 */
1158
1159 /* here we perform the actual fit: we pass now the median points via:
1160 * vec_wave,pvec_flux,nsamples params
1161 *
1162 * To fit each point in each fit window region we would instead pass
1163 * pwav,pflux,nsel params
1164 *
1165 */
1166
1167 double* pvec_wave=NULL;
1168 double* pvec_flux=NULL;
1169 pvec_wave=cpl_vector_get_data(vec_wave);
1170 pvec_flux=cpl_vector_get_data(vec_flux);
1171 //cpl_vector_save(vec_wave,"papero.fits",CPL_BPP_IEEE_FLOAT,NULL,CPL_IO_DEFAULT);
1172 //xsh_msg("Max wave=%g",cpl_vector_get_max(vec_wave));
1173
1174 //sfit = xsh_bspline_fit_data2(pvec_wave, pvec_flux, nsamples, phigh,instrument, 1);
1175
1176
1177 check(sfit = xsh_bspline_interpolate_data_at_pos(pvec_wave, pvec_flux, nsamples+2,
1178 pvec_wave,nsamples+2));
1179
1180
1181 //xsh_msg("Max wave=%g",cpl_vector_get_max(vec_wave));
1182
1183 //xsh_free_vector(&vflux_med);
1184 //cpl_vector_unwrap(vflux);
1185 int i=0;
1186 /* this instruction is to fill the fit column in table tab_stack with result of the fit
1187 for(i=0;i<nsel;i++) {
1188 pfit[i]=sfit[i];
1189 xsh_msg("sfit=%g",sfit[i]);
1190 }
1191
1192
1193 sprintf(fname,"ratio_fit.fits");
1194 check(xsh_table_save(tab_stack,NULL,NULL,fname,ext));
1195 */
1196 xsh_free_table(&tab_stack);
1197
1198 cpl_table* tab_med = cpl_table_new(nsamples+2);
1199
1200 check(cpl_table_wrap_double(tab_med,pvec_wave,"wave"));
1201 check(cpl_table_wrap_double(tab_med,pvec_flux,"ratio"));
1202 check(cpl_table_wrap_double(tab_med,sfit,"fit"));
1203 sprintf(fname,"ratio_med.fits");
1204 //check(xsh_table_save(tab_med,NULL,NULL,fname,ext));
1205 //xsh_msg("Max wave=%g",cpl_vector_get_max(vec_wave));
1206 /*
1207 check(
1208 tab_corr=xsh_table_downsample_table(tab_med,"wave","fit",table_rs,"wave","fit"));
1209 */
1210 /*
1211 int nrow=cpl_table_get_nrow(tab_corr);
1212 */
1213 double* wave=NULL;
1214
1215 tab_corr=cpl_table_duplicate(table_rs);
1216 int nrow=cpl_table_get_nrow(tab_corr);
1217 cpl_table_new_column(tab_corr,"fit",CPL_TYPE_DOUBLE);
1218 cpl_table_fill_column_window_double(tab_corr, "fit", 0, nrow, 0);
1219
1220 check(pfit=cpl_table_get_data_double(tab_corr,"fit"));
1221 check(wave=cpl_table_get_data_double(tab_corr,"wave"));
1222 //xsh_msg("Max wave=%g",cpl_vector_get_max(vec_wave));
1223 cpl_free(sfit);
1224 check(sfit = xsh_bspline_interpolate_data_at_pos(pvec_wave, pvec_flux,
1225 nsamples+2,wave,nrow));
1226 for (i = 0; i < nrow; i ++) {
1227 pfit[i]=sfit[i];
1228 }
1229 //check(cpl_table_wrap_double(tab_corr,sfit,"fit"));
1230
1231 //check(xsh_table_save(tab_corr,NULL,NULL,"ratio_fit_continuum.fits",ext));
1232
1233 //check(xsh_table_save(table_rs,NULL,NULL,"table_rs_check.fits",ext));
1234
1235 /*
1236 9) divide telluric corrected spectrum by continuum fit
1237 */
1238 check(cpl_table_duplicate_column(tab_corr, "correction", table_rs, "ratio"));
1239 check(cpl_table_divide_columns(tab_corr, "correction", "fit"));
1240 //sprintf(fname,"spectrum_obs_corrected.fits");
1241
1242
1243 cpl_table_unwrap(tab_med,"wave");
1244 cpl_table_unwrap(tab_med,"ratio");
1245 cpl_table_unwrap(tab_med,"fit");
1246 xsh_free_table(&tab_med);
1247 xsh_free_table(&table_rs);
1248 xsh_free_vector(&vec_wave);
1249 xsh_free_vector(&vec_flux);
1250
1251 cpl_free(sfit);
1252
1253
1254 /*
1255 10) find a way to identify the best correction from the residuals
1256 (I tried some statistical tests, but am not really convinced)
1257 */
1258 check(xsh_evaluate_tell_model(tab_corr,instrument,ext,&mean,&rms));
1259 qc_head=cpl_propertylist_new();
1260 cpl_propertylist_append_double(qc_head,XSH_QC_TELLCORR_RATAVG,mean);
1261 cpl_propertylist_append_double(qc_head,XSH_QC_TELLCORR_RATRMS,rms);
1262 //cpl_propertylist_dump(qc_head,stdout);
1263 //check(xsh_table_save(tab_corr,NULL,qc_head,fname,ext));
1264 xsh_msg("ext=%d mean=%g rms=%g", ext,mean, rms);
1265 xsh_free_propertylist(&qc_head);
1266 xsh_free_table(&ext_tab);
1267 xsh_free_table(&tab_corr);
1268 pmean[ext]=mean;
1269 prms[ext]=rms;
1270 pext[ext]=ext;
1271 xsh_free_table(&table_rs_div_mc);
1272 }
1273 double model_mean=0;
1274 double model_rms=0;
1275 int status=0;
1276 check(model_mean=cpl_table_get_column_min(tab_res,"mean"));
1277 //check(xsh_table_save(tab_res,NULL,NULL,"tab_res.fits",0));
1278 //xsh_msg("optimal model mean %g",model_mean);
1279 cpl_table_get_column_minpos(tab_res,"mean",model_idx);
1280 check(model_rms=cpl_table_get_double(tab_res,"rms",*model_idx,&status));
1281 xsh_msg("optimal model: mean=%g rms=%g",model_mean,model_rms);
1282 xsh_msg("optimal model id %d mean %g rms %g",
1283 (int)(*model_idx),model_mean,model_rms);
1284
1285 /* load finally optimal solution */
1286 //xsh_msg("fname=%s",fname);
1287 sprintf(fname,"spectrum_observed_divided_by_model_convolved.fits");
1289 ext=(int)(*model_idx);
1290 xsh_msg("fname=%s ext=%d",fname,ext);
1291 check(table_rs=cpl_table_load(fname,ext+1,0));
1292 //xsh_msg("table_rs nrow=%d",(int)cpl_table_get_nrow(table_rs));
1293 qc_head=cpl_propertylist_new();
1294 cpl_propertylist_append_double(qc_head,XSH_QC_TELLCORR_RATAVG,mean);
1295 cpl_propertylist_append_double(qc_head,XSH_QC_TELLCORR_RATRMS,rms);
1296 cpl_propertylist_append_int(qc_head,XSH_QC_TELLCORR_OPTEXTID,ext+1);
1297 //cpl_propertylist_dump(qc_head,stdout);
1298 //sprintf(fname,"model_optimal.fits");
1299 check(xsh_table_save(table_rs,NULL,qc_head,fname,0));
1300 //cpl_table_dump(table_rs,1,2,stdout);
1301
1302 //xsh_msg("table_rm nrow=%d",(int)cpl_table_get_nrow(table_rm));
1303 //cpl_table_dump(table_rm,1,2,stdout);
1304 //xsh_msg("ok1");
1305
1306 //check(xsh_table_save(table_rs,NULL,NULL,"check.fits",0));
1307 check(tab_corr=xsh_table_downsample_table(table_rs,"wave","ratio",table_rm,"wavelength","ratio"));
1308 //xsh_msg("tab_cor nrow=%d",(int)cpl_table_get_nrow(tab_corr));
1309 //cpl_table_dump(tab_corr,1,2,stdout);
1310 //check(xsh_table_save(tab_corr,NULL,NULL,"pippo.fits",0));
1311
1312 cleanup:
1313 xsh_free_propertylist(&qc_head);
1314 xsh_free_table(&table_rs);
1315 xsh_free_table(&tab_res);
1316 xsh_free_table(&table_s);
1317 xsh_free_table(&table_rm);
1318
1319 return tab_corr;
1320
1321}
1322
1323
static double sigma
static xsh_instrument * instrument
int xsh_spectrum_get_size_lambda(xsh_spectrum *s)
Get lambda axis size of spectrum.
double xsh_spectrum_get_lambda_min(xsh_spectrum *s)
Get minimum lambda of spectrum.
double * xsh_spectrum_get_errs(xsh_spectrum *s)
Get errs of spectrum.
double xsh_spectrum_get_lambda_step(xsh_spectrum *s)
Get bin in lambda of spectrum.
int * xsh_spectrum_get_qual(xsh_spectrum *s)
Get qual of spectrum.
double * xsh_spectrum_get_flux(xsh_spectrum *s)
Get flux of spectrum.
double xsh_spectrum_get_lambda_max(xsh_spectrum *s)
Get maximum lambda of spectrum.
#define check(COMMAND)
Definition: xsh_error.h:71
XSH_ARM xsh_instrument_get_arm(xsh_instrument *i)
Get an arm on instrument structure.
int size
static SimAnneal s
Definition: xsh_model_sa.c:99
#define xsh_msg(...)
Print a message on info level.
Definition: xsh_msg.h:121
static cpl_error_code xsh_align_model_to_spectrum(cpl_table *table_me, cpl_table *table_se, const double wlogstp, const int norm_xcorr, const int ext, cpl_table **table_mm)
static cpl_error_code xsh_extract_points_to_fit(cpl_table *table, const char *cwav, const char *cratio, const int nsamples, xsh_instrument *instrument, cpl_vector **vec_wave, cpl_vector **vec_flux)
cpl_table * xsh_table_resample_uniform(cpl_table *tinp, const char *cwinp, const char *cfinp, const double wstp)
static double xsh_resample_double(double wout, double *pw1, double *pf1, double wmin, double wmax, int size_obs)
static cpl_table * xsh_spectrum_to_table(xsh_spectrum *s)
cpl_error_code xsh_correl_spectra(double *flux_s, double *flux_m, const int size, const int hsearch, const double wlogstp, const int norm_xcorr, const double range, const int ext, XSH_GAUSSIAN_FIT *gfit)
cpl_table * xsh_telluric_model_eval(cpl_frame *frame_m, xsh_spectrum *s, xsh_instrument *instrument, cpl_size *model_idx)
static cpl_table * xsh_table_select_range(cpl_table *table_in, const char *col, const double wmin, const double wmax)
static cpl_table * xsh_table_resample_table(cpl_table *tinp, const char *cwinp, const char *cfinp, cpl_table *tref, const char *cwref, const char *cfref)
static cpl_error_code xsh_evaluate_tell_model(cpl_table *corr, xsh_instrument *instrument, const int ext, double *mean, double *rms)
static cpl_error_code xsh_get_xcorrel_peak(cpl_vector *wcorr, cpl_vector *fcorr, XSH_GAUSSIAN_FIT *gfit, const double range, const int ext)
static cpl_table * xsh_extract_ranges_to_fit(cpl_table *table_in, const char *cwav, xsh_instrument *instrument)
static cpl_table * xsh_table_downsample_table(cpl_table *tinp, const char *cwinp, const char *cfinp, cpl_table *tref, const char *cwref, const char *cfref)
void xsh_free_vector(cpl_vector **v)
Deallocate a vector and set the pointer to NULL.
Definition: xsh_utils.c:2284
double * xsh_function1d_xcorrelate(double *line_i, int width_i, double *line_t, int width_t, int half_search, int normalise, double *xcorr_max, double *delta)
Definition: xsh_utils.c:7101
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
double * xsh_bspline_interpolate_data_at_pos(double *w_data, double *f_data, const int n_data, double *w_pos, const int n_pos)
@ XSH_ARM_VIS
int m
Definition: xsh_detmon_lg.c:91
cpl_error_code xsh_table_save(cpl_table *t, cpl_propertylist *ph, cpl_propertylist *xh, const char *fname, const int ext)
Definition: xsh_dfs.c:5171
HIGH_ABS_REGION * xsh_fill_tell_compute_resid_regions(xsh_instrument *instrument, cpl_frame *high_abs_frame)
HIGH_ABS_REGION * xsh_fill_tell_fit_regions(xsh_instrument *instrument, cpl_frame *high_abs_frame)
#define XSH_QC_TELLCORR_RATRMS
Definition: xsh_pfits_qc.h:64
#define XSH_QC_TELLCORR_OPTEXTID
Definition: xsh_pfits_qc.h:65
#define XSH_QC_TELLCORR_RATAVG
Definition: xsh_pfits_qc.h:63