X-shooter Pipeline Reference Manual 3.8.15
xsh_utils_efficiency.c
Go to the documentation of this file.
1/*
2 * This file is part of the ESO X-Shooter Pipeline
3 * Copyright (C) 2004-2009 European Southern Observatory
4 *
5 * This program 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-04-22 08:41:32 $
22 * $Revision: 1.59 $
23 *
24 */
25
26#ifdef HAVE_CONFIG_H
27#include <config.h>
28#endif
29#include <string.h>
30#include <math.h>
31#include <cpl.h>
32/* irplib */
33#include <irplib_utils.h> /* portable isinf check */
34
35
36#include <xsh_pfits.h>
37#include <xsh_msg.h>
38#include <xsh_dfs.h>
39#include <xsh_error.h>
40#include <xsh_utils.h>
41#include <xsh_utils_table.h>
42#include <xsh_utils_wrappers.h>
45#include <xsh_star_index.h>
46#include <xsh_data_atmos_ext.h>
47
48#define PRO_STD_STAR_SPECTRA "STD_STAR_SPECTRA"
49static const char COL_NAME_HIGH_ABS[] = "HIGH_ABS";
50static const char COL_NAME_WAVELENGTH[] = "WAVELENGTH";
51//static const char COL_NAME_WAVELENGTH_C[] = "WAVELENGTH";
52static const char COL_NAME_WAVE_ATMDISP[] = "LAMBDA";
54
55static const char COL_NAME_REF[] = "REF";
56static const char COL_NAME_COR[] = "COR";
57static const char COL_NAME_SRC_COR[] = "SRC_COR";
58
59static const char COL_NAME_WAVE_OBJ[] = "WAVELENGTH";
60static const char COL_NAME_INT_OBJ[] = "INT_OBJ";//"counts_tot";
61static const char COL_NAME_ORD_OBJ[] = "ORD";//"counts_tot";
62static const char COL_NAME_WAVE_REF[] = "LAMBDA";
63static const char COL_NAME_FLUX_REF[] = "FLUX";
64static const char COL_NAME_BINWIDTH_REF[] = "BIN_WIDTH";
65static const char COL_NAME_EPHOT[] = "EPHOT";
66static const char COL_NAME_EXT[] = "EXT";
67static const char COL_NAME_SRC_EFF[] = "EFF";
68//static const char PAR_NAME_DIT[] = "ESO DET DIT";
69//static const double UVES_flux_factor =1.e16;
70
72
73static int
74xsh_column_to_double(cpl_table* ptable, const char* column);
75
76
77static double*
78xsh_create_column_double(cpl_table* tbl, const char* col_name, int nrow);
79
80
81
82static cpl_error_code
83xsh_get_std_obs_values(cpl_propertylist* plist,
84 double* exptime,
85 double* airmass,
86 double* dRA,
87 double* dDEC);
88
89/*--------------------------------------------------------------------------*/
92/*---------------------------------------------------------------------------*/
104/*---------------------------------------------------------------------------*/
105
106void
107xsh_load_ref_table(cpl_frameset* frames,
108 double dRA,
109 double dDEC,
110 double EPSILON,
112 cpl_table** pptable)
113{
114 const char* name = NULL;
115 cpl_frame* frm_ref = NULL;
116
117 /* REF STD frame */
119 if (!frm_ref)
120 {
121 xsh_msg("REF frame is not found, trying to get REF from the catalog");
122 /* REF STD catalog frame */
123 // get from catalog
125 if (frm_ref)
126 {
127 check(name=cpl_frame_get_filename(frm_ref));
128 if (name)
129 {
130 star_index* pstarindex = star_index_load(name);
131 if (pstarindex)
132 {
133 const char* star_name = NULL;
134 xsh_msg("Searching std RA[%f] DEC[%f] with tolerance[%f] in star catalog", dRA, dDEC, EPSILON);
135 *pptable = star_index_get(pstarindex, dRA, dDEC, EPSILON, EPSILON, &star_name);
136
137 if (*pptable && star_name)
138 {
139 xsh_msg("Found STD star: %s", star_name);
140 }
141 else
142 {
143 xsh_msg("ERROR - REF star %s could not be found in the catalog",star_name);
144 }
145 }
146 else
147 {
148 xsh_msg("ERROR - could not load the catalog");
149 }
150 }
151 }
152 }
153 else
154 {
155 xsh_msg("REF frame is found");
156 check(name=cpl_frame_get_filename(frm_ref));
157 check(*pptable=cpl_table_load(name,1,0));
158 }
159 return;
160 cleanup:
161 return;
162}
163
164
165cpl_error_code
167
168 switch (arm) {
169 case XSH_ARM_AGC:
170 case XSH_ARM_UVB:
171 w->wguess = 486.1322;
172 w->range_wmin=450;
173 w->range_wmax=550;
174 w->ipol_wmin_max=470;
175 w->ipol_wmax_min=500;
176 w->ipol_hbox=0.5;
177
178 break;
179 case XSH_ARM_VIS:
180 w->wguess= 656.2793;
181 w->range_wmin=640;
182 w->range_wmax=670;
183 w->ipol_wmin_max=653;
184 w->ipol_wmax_min=661;
185 w->ipol_hbox=1.5;
186
187 break;
188 case XSH_ARM_NIR:
189 w->wguess= 1094.12;
190 w->range_wmin=1000;
191 w->range_wmax=1200;
192 w->ipol_wmin_max=1090;
193 w->ipol_wmax_min=1100;
194 w->ipol_hbox=1.5;
195
196 break;
198
199 w->wguess= 999;
200 w->range_wmin=999;
201 w->range_wmax=999;
202 w->ipol_wmin_max=999;
203 w->ipol_wmax_min=999;
204 w->ipol_hbox=1.5;
205 break;
206 }
207
208 return cpl_error_get_code();
209}
210
212
213 xsh_rv_ref_wave_param * self = cpl_malloc(sizeof(xsh_rv_ref_wave_param));
214 self->wguess=999.; /* Reference line wavelength position */
215 self->range_wmin=999.; /* minimum of wavelength box for line fit */
216 self->range_wmax=999.; /* maximum of wavelength box for line fit */
217 self->ipol_wmin_max=999.; /* maximum lower sub-range wavelength value used to fit line slope */
218 self->ipol_wmax_min=999.; /* minimum upper sub-range wavelength value used to fit line slope */
219
220 return self;
221}
222
224
225 cpl_free(p);
226
227 return ;
228}
229
230/*---------------------------------------------------------------------------*/
241/*---------------------------------------------------------------------------*/
242
243cpl_error_code xsh_parse_catalog_std_stars(cpl_frame* cat, double dRA,
244 double dDEC, double EPSILON, cpl_table** pptable,xsh_std_star_id* std_star_id) {
245 const char* name = NULL;
246 XSH_ASSURE_NOT_NULL_MSG(cat, "Provide input catalog");
247 if (cat) {
248 check(name=cpl_frame_get_filename(cat));
249 if (name) {
250 star_index* pstarindex = star_index_load(name);
251 if (pstarindex) {
252 const char* star_name = NULL;
253 xsh_msg(
254 "Searching std RA[%f] DEC[%f] with tolerance[%f] in star catalog", dRA, dDEC, EPSILON);
255 *pptable = star_index_get(pstarindex, dRA, dDEC, EPSILON, EPSILON,
256 &star_name);
257 if(star_name != NULL) {
258 if ( strcmp(star_name,"GD71") == 0 ) *std_star_id = XSH_GD71;
259 else if ( strcmp(star_name,"Feige110") == 0 ) *std_star_id = XSH_Feige110;
260 else if ( strcmp(star_name,"GD153") == 0 ) *std_star_id = XSH_GD153;
261 else if ( strcmp(star_name,"LTT3218") == 0 ) *std_star_id = XSH_LTT3218;
262 else if ( strcmp(star_name,"LTT7987") == 0 ) *std_star_id = XSH_LTT7987;
263 else if ( strcmp(star_name,"BD17") == 0 ) *std_star_id = XSH_BD17;
264 else if ( strcmp(star_name,"EG274") == 0 ) *std_star_id = XSH_EG274;
265 }
266 xsh_msg("star index=%d",*std_star_id);
267 if (*pptable && star_name != NULL) {
268 xsh_msg("Found STD star: %s", star_name);
269
270 } else {
271 xsh_msg(
272 "ERROR - REF star %s could not be found in the catalog", star_name);
273 }
274 } else {
275 xsh_msg("ERROR - could not load the catalog");
276 }
277 star_index_delete(pstarindex);
278 }
279 }
280 cleanup: return cpl_error_get_code();
281}
282
283
284
285
286
287static double*
288xsh_create_column_double(cpl_table* tbl, const char* col_name, int nrow)
289{
290 double* retval = 0;
291 check(cpl_table_new_column(tbl, col_name, CPL_TYPE_DOUBLE));
292 check(cpl_table_fill_column_window_double(tbl, col_name, 0, nrow, -1));
293 check(retval = cpl_table_get_data_double(tbl,col_name));
294 return retval;
295 cleanup:
296 return retval;
297}
298
299
300/*---------------------------------------------------------------------------*/
310/*---------------------------------------------------------------------------*/
311
312cpl_error_code
313xsh_get_std_obs_values(cpl_propertylist* plist,
314 double* exptime,
315 double* airmass,
316 double* dRA,
317 double* dDEC)
318{
319
320 *airmass=xsh_pfits_get_airm_mean(plist) ;
321 *dDEC=xsh_pfits_get_dec(plist);
322 *dRA=xsh_pfits_get_ra(plist);
324
325 return cpl_error_get_code();
326
327}
328
329cpl_error_code
331{
332 int nrow=cpl_table_get_nrow(*eff);
333 double* pwav=NULL;
334 int* phigh_abs=NULL;
335 int k=0;
336 int i=0;
337 cpl_table_new_column(*eff,COL_NAME_HIGH_ABS,CPL_TYPE_INT);
338 cpl_table_fill_column_window_int(*eff,COL_NAME_HIGH_ABS,0,nrow,0);
339 pwav=cpl_table_get_data_double(*eff,COL_NAME_WAVELENGTH);
340 phigh_abs=cpl_table_get_data_int(*eff,COL_NAME_HIGH_ABS);
341 if ( phigh != NULL ) {
342 for( k = 0 ; phigh->lambda_min != 0. ; k++, phigh++ ) {
343 for(i=0;i<nrow;i++) {
344 if(pwav[i]>=phigh->lambda_min && pwav[i]<=phigh->lambda_max) {
345 phigh_abs[i]=1;
346 }
347 }
348 }
349 }
350
351 return cpl_error_get_code();
352
353}
354
355/*---------------------------------------------------------------------------*/
373/*---------------------------------------------------------------------------*/
374
375cpl_frame*
377 cpl_frameset * frames,
378 double dGain,
379 double dEpsilon,
380 double aimprim,
381 xsh_instrument* inst,
382 const char* col_name_atm_wave,
383 const char* col_name_atm_abs,
384 const char* col_name_ref_wave,
385 const char* col_name_ref_flux,
386 const char* col_name_ref_bin,
387 const char* col_name_obj_wave,
388 const char* col_name_obj_flux
389 )
390{
391 cpl_frame* frm_sci = NULL;
392 cpl_frame* frm_atmext = NULL;
393 cpl_table* tbl_obj_spectrum = NULL; // input table with spectrum
394 cpl_table* tbl_atmext = NULL;
395 double exptime = 600;
396
397 cpl_propertylist* plist = NULL;
398 cpl_table* tbl_ref = NULL;
399 cpl_frame* result=NULL;
400
401 const char* name=NULL;
402 double dRA = 0;
403 double dDEC = 0;
404 char prod_name[256];
405 char prod_tag[256];
406 double airmass=0;
407 const double mk2AA=1E4; /* mkm/AA */
408 int nclip=0;
409 int ntot=0;
410 cpl_table* tbl_result = NULL;
411
412
413 // read all input data and call the internal function
414 // read the input tables
415 // 1. observation
417 check(name=cpl_frame_get_filename(frm_sci));
418 //xsh_msg("name=%s",name);
419 check(tbl_obj_spectrum=cpl_table_load(name,1,0));
420 check(plist=cpl_propertylist_load(name,0));
421
422 xsh_get_std_obs_values(plist,&exptime,&airmass,&dRA,&dDEC);
423
424 // 2. reference table
425 xsh_load_ref_table(frames, dRA, dDEC, dEpsilon, inst, &tbl_ref);
426 if (tbl_ref)
427 {
428 // 3. atmext
429 check(frm_atmext=cpl_frameset_find(frames,FRM_EXTCOEFF_TAB));
430 check(name=cpl_frame_get_filename(frm_atmext));
431 check(tbl_atmext=cpl_table_load(name,1,0));
432
434 tbl_obj_spectrum,
435 tbl_atmext,
436 tbl_ref,
437 exptime,
438 airmass,
439 aimprim,
440 dGain,
441 1,
442 mk2AA,//valid only for SINFONI
443 col_name_atm_wave,
444 col_name_atm_abs,
445 col_name_ref_wave,
446 col_name_ref_flux,
447 col_name_ref_bin,
448 col_name_obj_wave,
449 col_name_obj_flux,
450 &ntot,&nclip);
451 if (tbl_result)
452 {
453 HIGH_ABS_REGION * phigh=NULL;
454 check(xsh_efficiency_add_high_abs_regions(&tbl_result,phigh));
455 sprintf(prod_tag,"EFFICIENCY_%s",xsh_instrument_arm_tostring(inst));
456 sprintf(prod_name,"%s.fits",prod_tag);
457
458 result=xsh_frame_product(prod_name,prod_tag,CPL_FRAME_TYPE_TABLE,
459 CPL_FRAME_GROUP_CALIB,CPL_FRAME_LEVEL_FINAL);
460 cpl_table_save(tbl_result, plist, NULL,prod_name, CPL_IO_DEFAULT);
461
462 xsh_free_table(&tbl_result);
463 }
464
465 }
466
467 cleanup:
468 xsh_free_propertylist(&plist);
469 xsh_free_table(&tbl_atmext);
470 xsh_free_table(&tbl_obj_spectrum);
471 xsh_free_table(&tbl_ref);
472 return result;
473}
474
475static int
476xsh_column_to_double(cpl_table* ptable, const char* column)
477{
478 const char* TEMP = "_temp_";
479 check(cpl_table_duplicate_column(ptable, TEMP, ptable, column));
480 check(cpl_table_erase_column(ptable, column));
481 check(cpl_table_cast_column(ptable, TEMP, column, CPL_TYPE_DOUBLE));
482 check(cpl_table_erase_column(ptable, TEMP ));
483 return 0;
484 cleanup:
485 xsh_msg(" error column to double [%s]", column);
486 return -1;
487}
488
489
490/*---------------------------------------------------------------------------*/
512/*---------------------------------------------------------------------------*/
513
514cpl_table*
515xsh_utils_efficiency_internal(cpl_table* tbl_obj_spectrum,
516 cpl_table* tbl_atmext,
517 cpl_table* tbl_ref,
518 double exptime,
519 double airmass,
520 double aimprim,
521 double gain,
522 int biny,
523 double src2ref_wave_sampling,
524 const char* col_name_atm_wave,
525 const char* col_name_atm_abs,
526 const char* col_name_ref_wave,
527 const char* col_name_ref_flux,
528 const char* col_name_ref_bin,
529 const char* col_name_obj_wave,
530 const char* col_name_obj_flux,
531 int* ntot,
532 int* nclip
533 )
534{
535
536 const double TEL_AREA = 51.2e4;//49.8340828812e4;// /* collecting area in cm2 */
537 double cdelta1 = 0;
538 int i = 0;
539 cpl_table* tbl_sel = NULL;
540
541 /* structure of the input table
542 * col type description
543 * wavelength double
544 * flux double
545 * */
546 cpl_table* tbl_result = NULL; // output table
547
548 /*
549 * structure of the output table
550 * col type description
551 * wavelength double
552 * EFF double efficiency in range 0-1
553 * */
554
555 double* pref = NULL;
556 double* pext = NULL;
557 double* pcor = NULL;
558 double* peph = NULL;
559
560 double* pw = NULL; // wavelength of the input table
561 int nrow = 0;
562 double ref_bin_size=0;
563 double um2nm=0.001;
564 double eff_med=0.;
565 double eff_rms=0.;
566 double eff_thresh=0.;
567
568 double kappa=5;
569 double nm2AA=10.;
570 int nsel=0;
571 double eff=0;
572 nrow = cpl_table_get_nrow(tbl_obj_spectrum);
573 xsh_msg_dbg_medium("Starting efficiency calculation: exptime[%f] airmass[%f] nrow[%d]",
574 exptime, airmass, nrow);
575
576 //cpl_table_dump(tbl_obj_spectrum,0,3,stdout);
577 //xsh_msg("col wave=%s col_flux=%s",col_name_obj_wave,col_name_obj_flux);
578
579 /* convert to double */
580 xsh_column_to_double(tbl_obj_spectrum,col_name_obj_wave);
581 xsh_column_to_double(tbl_obj_spectrum,col_name_obj_flux);
582
583 check(xsh_column_to_double(tbl_atmext,col_name_atm_wave ));
584 check(xsh_column_to_double(tbl_atmext,col_name_atm_abs ));
585 check(xsh_column_to_double(tbl_ref,col_name_ref_wave ));
586 check(xsh_column_to_double(tbl_ref,col_name_ref_flux ));
587 check(xsh_column_to_double(tbl_ref,col_name_ref_bin ));
588
589 /* as the reference spectra are erg cm-2 s-1 AA-1 and the wave scale
590 have been converted from AA to nm to be uniform to XSH scale
591 we apply a corresponding factor to the flux */
592 //check(cpl_table_multiply_scalar(tbl_ref,col_name_ref_flux,nm2AA));
593
594 /* get bin size of ref STD star spectrum */
595 ref_bin_size=cpl_table_get_double(tbl_ref,col_name_ref_bin,0,NULL);
596 xsh_msg_dbg_medium("ref_bin_size[nm]=%g",ref_bin_size);
597
598 /*
599 xsh_msg("Unit conversion factor src/ref STD spectrum=%g",
600 src2ref_wave_sampling);
601 */
602
603 /* convert obj spectrum wave unit to the same of the reference one
604 This is required to be able to obtain the interpolated value of
605 the extinction and reference std star spectrum at the sampled
606 object wavelength when we compute the correction term due to
607 atmospheric extinction
608 */
609 check(cpl_table_multiply_scalar(tbl_obj_spectrum,col_name_obj_wave,
610 src2ref_wave_sampling));
611
612
613 /* Computes the atmospheric distorsion correction:
614 COL_NAME_REF: reference std star flux value
615 COL_NAME_EXT: atm extinction
616 COL_NAME_EPHOT: photon energy
617 COL_NAME_COR: corrective factor (due to atmospheric extinction)
618 */
619 xsh_msg_dbg_medium("object 2 src std wave factor %g",src2ref_wave_sampling);
620 check(pw=cpl_table_get_data_double(tbl_obj_spectrum,col_name_obj_wave));
621
622 // prepare columns for the output data
623 check(tbl_result=cpl_table_new(nrow));
624 check(pref=xsh_create_column_double(tbl_result, COL_NAME_REF, nrow));
625 check(pext=xsh_create_column_double(tbl_result, COL_NAME_EXT, nrow));
626 check(pcor=xsh_create_column_double(tbl_result, COL_NAME_COR, nrow));
627 check(peph=xsh_create_column_double(tbl_result, COL_NAME_EPHOT, nrow));
628 xsh_msg_dbg_medium("wave range: [%g,%g] nm",pw[0],pw[nrow-1]);
629 xsh_msg_dbg_medium("src_bin_size[nm]=%g",pw[1] - pw[0]);
630
631 cdelta1 = (pw[1] - pw[0]) / src2ref_wave_sampling ; /* we want the delta in original units. As we rescaled to AA we need to correct for this */
632 xsh_msg_dbg_medium("nrow=%d cdelta1=%g",nrow,cdelta1);
633 for (i = 0; i < nrow; i++)
634 {
635 check(pext[i] = xsh_table_interpolate(tbl_atmext, pw[i],col_name_atm_wave, col_name_atm_abs));
636 check(pref[i] = xsh_table_interpolate(tbl_ref, pw[i], col_name_ref_wave,col_name_ref_flux));
637 pcor[i] = pow(10,(-0.4*pext[i] * (aimprim - airmass)));
638 eff = 1.e7*1.986e-19/(pw[i]*um2nm);
639 if(!isnan(eff)) {
640 peph[i]=eff;
641 } else {
642 /*
643 xsh_msg("pw=%g pext=%g pref=%g pcor=%g pef=%g",
644 pw[i],pext[i],pref[i],pcor[i],peph[i]);
645 */
646 cpl_table_set_invalid(tbl_result,COL_NAME_EPHOT,i);
647 }
648 /* ph energy: 1.986*10^-19(J.ph^-1)/lam(um) ==>
649 in as pw is expressed in nm units we need to multiply by um2nm: 10-3
650 */
651 /*
652 if(i< 2) {
653 xsh_msg("pw[i]=%g,pcor=%g peph=%g",pw[i],pcor[i],peph[i]);
654 }
655 */
656 }
657 cpl_table_erase_invalid(tbl_result);
658
659 /*
660 xsh_msg("atm: %s, %s ref: %s, %s obj: %s, %s",
661 col_name_atm_wave,col_name_atm_abs,
662 col_name_ref_wave,col_name_ref_flux,
663 col_name_obj_wave,col_name_obj_flux);
664 */
665 /* add in result table also ORDER, OBJ(flux), WAVELENGTH columns */
666 check(cpl_table_duplicate_column(tbl_result,"ORDER",
667 tbl_obj_spectrum,COL_NAME_ORD_OBJ ));
668 check(cpl_table_duplicate_column(tbl_result,COL_NAME_SRC_COR,
669 tbl_obj_spectrum, col_name_obj_flux));
670 check(cpl_table_duplicate_column(tbl_result,col_name_obj_wave,
671 tbl_obj_spectrum,col_name_obj_wave));
672
673 /* correct for atmospheric extintion:
674 bring the observed STD star out of the atmosphere */
675 check(cpl_table_multiply_columns(tbl_result,COL_NAME_SRC_COR,COL_NAME_COR));
676
677 /* correct object flux by binning:
678 As the reference std star spectra flux values are in units of
679 erg/cm2/s/A we have to obtain the bin size in Angstrom units.
680 Because the sampling step is in nm units, we need to apply a nm2AA factor:
681 cdelt1[src_sampling]*src2ref_wave_sampling*nm2AA */
682 cpl_table_divide_scalar(tbl_result, COL_NAME_SRC_COR, src2ref_wave_sampling);
683 cpl_table_divide_scalar(tbl_result, COL_NAME_SRC_COR, cdelta1);
684 cpl_table_divide_scalar(tbl_result, COL_NAME_SRC_COR, nm2AA);
685 /* correct for spatial bin size: Why?? */
686 cpl_table_divide_scalar(tbl_result,COL_NAME_SRC_COR,biny);
687
688
689 /*correct ref std star flux by binning:
690 divides by the bin size in ref_wave_sampling (usually AA) */
691 check(cpl_table_divide_scalar(tbl_result,COL_NAME_REF,ref_bin_size));
692
693 /* create efficiency column, initialized with the intensity value (ADU)
694 corrected for extinction and sampling step */
695 check(cpl_table_duplicate_column(tbl_result,COL_NAME_SRC_EFF,
696 tbl_result,COL_NAME_SRC_COR));
697
698 /* We need now to convert ADU to erg/cm2/s :
699 correct for detector gain, exposure time, telescope area: we are now
700 in units of electrons/s/cm2 */
701 check(cpl_table_multiply_scalar(tbl_result,COL_NAME_SRC_EFF,
702 gain / (exptime * TEL_AREA)));
703
704 /* To pass from electron to ergs we correct for photon energy */
705 check(cpl_table_multiply_columns(tbl_result,COL_NAME_SRC_EFF,
707
708 /* our observed STD star value is now finally in the same units
709 as the one of the reference STD star: the efficiency is the ratio
710 observed_flux/reference_flux */
711
712 check(cpl_table_divide_columns(tbl_result,COL_NAME_SRC_EFF,COL_NAME_REF));
713 /* apply factor UVES_flux_factor (1.e16) as reference catalog has fluxes
714 in units of UVES_flux_factor (1e-16) */
715 //check(cpl_table_multiply_scalar(tbl_result,COL_NAME_SRC_EFF,UVES_flux_factor));
716 //check(cpl_table_save(tbl_result,NULL,NULL,"pippo.fits", CPL_IO_DEFAULT));
717 /* To have cleaner plots we clean from outliers */
718 eff_med=cpl_table_get_column_median(tbl_result,COL_NAME_SRC_EFF);
719 eff_rms=cpl_table_get_column_stdev(tbl_result,COL_NAME_SRC_EFF);
720
721 eff_thresh=(eff_med+kappa*eff_rms<10.) ? eff_med+kappa*eff_rms:10.;
722 if(irplib_isinf(eff_thresh)) eff_thresh=10.;
723/*
724
725 cpl_table_and_selected_double(tbl_result,COL_NAME_SRC_EFF,
726 CPL_GREATER_THAN,1.e-5);
727*/
728/*
729 xsh_msg("1eff_med=%g eff_rms=%g thresh=%g",
730 eff_med,eff_rms,eff_thresh);
731*/
732 *ntot=cpl_table_get_nrow(tbl_result);
733 cpl_table_and_selected_double(tbl_result,COL_NAME_SRC_EFF,
734 CPL_GREATER_THAN,1.e-5);
735
736 cpl_table_and_selected_double(tbl_result,COL_NAME_SRC_EFF,
737 CPL_LESS_THAN,eff_thresh);
738
739 eff_med=cpl_table_get_column_median(tbl_result,COL_NAME_SRC_EFF);
740 eff_rms=cpl_table_get_column_stdev(tbl_result,COL_NAME_SRC_EFF);
741
742 eff_thresh=(eff_med+kappa*eff_rms<10.) ? eff_med+kappa*eff_rms:10.;
743 if(irplib_isinf(eff_thresh)) eff_thresh=10.;
744/*
745 xsh_msg("2eff_med=%g eff_rms=%g thresh=%g",
746 eff_med,eff_rms,eff_thresh);
747*/
748 cpl_table_and_selected_double(tbl_result,COL_NAME_SRC_EFF,
749 CPL_LESS_THAN,eff_thresh);
750
751 nsel=cpl_table_and_selected_double(tbl_result,COL_NAME_SRC_EFF,
752 CPL_LESS_THAN,1);
753
754 *nclip=(*ntot)-nsel;
755 tbl_sel=cpl_table_extract_selected(tbl_result);
756
757
758 cleanup:
759 xsh_free_table(&tbl_result);
760 return tbl_sel;
761}
762
771void
773 double* ra,
774 double* dec,
775 double* airmass)
776{
777
778 const char* name_sci=NULL;
779 cpl_propertylist* plist=NULL;
780
781 name_sci=cpl_frame_get_filename(frm_sci);
782 check(plist=cpl_propertylist_load(name_sci,0));
783 *ra=xsh_pfits_get_ra(plist);
784 *dec=xsh_pfits_get_dec(plist);
785 *airmass=xsh_pfits_get_airm_mean(plist);
786
787cleanup:
788
789 xsh_free_propertylist(&plist);
790 return;
791}
792
793
794// TO BE FIXED (gain is wrong in NIR)
795static void
798 double* gain,
799 double* airmass,
800 double* exptime,
801 int* naxis1,
802 int* biny)
803{
804
805 const char* name_sci=NULL;
806 cpl_propertylist* plist=NULL;
807
808 name_sci=cpl_frame_get_filename(frm_sci);
809 check(plist=cpl_propertylist_load(name_sci,0));
810 check(*naxis1=xsh_pfits_get_naxis1(plist));
811 *airmass=xsh_pfits_get_airm_mean(plist);
812
814 *gain=2.12;
815 *biny=1;
817 } else {
818 check(*gain=xsh_pfits_get_conad(plist));
819 check(*biny=xsh_pfits_get_biny(plist));
821 }
822
823cleanup:
824 xsh_free_propertylist(&plist);
825
826 return;
827}
828
834cpl_frame*
836 cpl_frame* frm_sci)
837{
838 cpl_frame* result=NULL;
839 double dRA=0;
840 double dDEC=0;
841 //double nm2AA=10.;
842 double nm2nm=1.;
843 cpl_table* tbl_ref=NULL;
844 char fname[256];
845 char ftag[256];
846 double airmass=0;
847 xsh_std_star_id std_star_id=0;
848
849 XSH_ASSURE_NOT_NULL_MSG(frm_sci,"Null input sci frame set!Exit");
850 XSH_ASSURE_NOT_NULL_MSG(frm_cat,"Null input std star cat frame set!Exit");
851 xsh_frame_sci_get_ra_dec_airmass(frm_sci,&dRA,&dDEC,&airmass);
852 check(xsh_parse_catalog_std_stars(frm_cat,dRA,dDEC,STAR_MATCH_DEPSILON,&tbl_ref,&std_star_id));
853
854 //cpl_table_divide_scalar(tbl_ref,COL_NAME_FLUX_REF,UVES_flux_factor);
855 cpl_table_divide_scalar(tbl_ref,COL_NAME_WAVE_REF,nm2nm);
856 cpl_table_multiply_scalar(tbl_ref,COL_NAME_FLUX_REF,nm2nm);
857 check(cpl_table_divide_columns(tbl_ref,COL_NAME_FLUX_REF,COL_NAME_BINWIDTH_REF));
858
859 sprintf(fname,"ref_std_star_spectrum.fits");
860 sprintf(ftag,"STD_STAR_FLUX");
861
862 check(cpl_table_save(tbl_ref,NULL,NULL,fname,CPL_IO_DEFAULT));
863 result=xsh_frame_product(fname,ftag,CPL_FRAME_TYPE_TABLE,
864 CPL_FRAME_GROUP_CALIB,CPL_FRAME_LEVEL_INTERMEDIATE);
865 cleanup:
866 return result;
867}
868
876cpl_frame*
877xsh_efficiency_compute(cpl_frame* frm_sci,
878 cpl_frame* frm_cat,
879 cpl_frame* frm_atmext,
880 cpl_frame* high_abs_win,
882
883{
884
885 cpl_image* ima_sci=NULL;
886 cpl_vector* vec_ord=NULL;
887 cpl_image* ima_obj=NULL;
888 cpl_table* obj_tab=NULL;
889
890 const char* name_sci=NULL;
891 const char* name_atm=NULL;
892 cpl_propertylist* plist=NULL;
893 cpl_propertylist* x_plist=NULL;
894
895
896 cpl_table* tbl_ord=NULL;
897 cpl_table* tot_eff=NULL;
898 cpl_table* tbl_eff=NULL;
899 cpl_table* tbl_ref=NULL;
900 cpl_table* tbl_atmext=NULL;
901
902 double * pobj=NULL;
903 double * pw=NULL;
904 double * pf=NULL;
905 int * po=NULL;
906
907 cpl_frame* frm_eff=NULL;
908
909 //char name_eff[256];
910
911 double crval1=0;
912 double cdelt1=0;
913 int naxis1=0;
914 int nrow=0;
915
916
917 double exptime=600;
918 cpl_vector* rec_profile=NULL;
919 int i=0;
920 double airmass=0;
921 double dRA=0;
922 double dDEC=0;
923 char key_name[40];
924 int next=0;
925 int nord=0;
926
927
928
929 int j=0;
930 double wav=0;
931 double gain=0;
932 int biny=1;
933 double aimprim=0;
934 int ord=0;
935 //double nm2AA=10.;
936 double nm2nm=1.;
937
938 char fname[256];
939 char tag[256];
940 double emax=0;
941 double emed=0;
942 int nclip=0;
943 int ntot=0;
944
945 int nclip_tot=0;
946 int neff_tot=0;
947 double fclip=0;
948 xsh_std_star_id std_star_id=0;
949 int vec_size=0;
950 XSH_ASSURE_NOT_NULL_MSG(frm_sci,"Null input sci frame set!Exit");
951 XSH_ASSURE_NOT_NULL_MSG(frm_cat,"Null input std star cat frame set!Exit");
952 XSH_ASSURE_NOT_NULL_MSG(frm_atmext,"Null input atmospheric ext frame set!Exit");
953
954
955 check(next = cpl_frame_get_nextensions(frm_sci));
956 //xsh_msg("next=%d",next);
957 nord=(next+1)/3;
958
959 xsh_frame_sci_get_ra_dec_airmass(frm_sci,&dRA,&dDEC,&airmass);
960 name_sci=cpl_frame_get_filename(frm_sci);
961 plist=cpl_propertylist_load(name_sci,0);
963 &gain,&airmass,&exptime,
964 &naxis1,&biny);
965
966
967
969 &tbl_ref,&std_star_id));
970
971
972
973 xsh_msg_dbg_medium("gain=%g airm=%g exptime=%g airmass=%g ra=%g dec=%g",
974 gain,airmass,exptime,airmass,dRA,dDEC);
975
976 xsh_msg_dbg_medium("name_sci=%s",name_sci);
977 nrow=naxis1*nord;
978 /* note that this table is overdimensioned: naxis1 is the size of the multi
979 extension image holding the 2D spectrum (usually order by order)
980 but orders may have different length==> different naxis1
981 */
982 obj_tab=cpl_table_new(nrow);
983
984
985 cpl_table_new_column(obj_tab,COL_NAME_ORD_OBJ,CPL_TYPE_INT);
986 cpl_table_new_column(obj_tab,COL_NAME_WAVE_OBJ,CPL_TYPE_DOUBLE);
987 cpl_table_new_column(obj_tab,COL_NAME_INT_OBJ,CPL_TYPE_DOUBLE);
988
989 check(cpl_table_fill_column_window_int(obj_tab,COL_NAME_ORD_OBJ,0,nrow,-1));
990 check(cpl_table_fill_column_window_double(obj_tab,COL_NAME_WAVE_OBJ,0,nrow,-1));
991 check(cpl_table_fill_column_window_double(obj_tab,COL_NAME_INT_OBJ,0,nrow,-1));
992
993
994
995 /* extract the object signal on a given slit
996 extract the sky signal on 2 given slits
997 computes the extracted sky per pixel on each sky window
998 average the two extracted sky
999 rescale the averaged sky to the object windows and
1000 subtract it from the object
1001 */
1002
1003 check(name_atm=cpl_frame_get_filename(frm_atmext));
1004 xsh_msg_dbg_medium("name_atm=%s",name_atm);
1005 check(tbl_atmext=cpl_table_load(name_atm,1,0));
1006
1007 if(!cpl_table_has_column(tbl_atmext,XSH_ATMOS_EXT_LIST_COLNAME_K)){
1008 xsh_msg_warning("You are using an obsolete atm extinction line table");
1009 cpl_table_duplicate_column(tbl_atmext,XSH_ATMOS_EXT_LIST_COLNAME_K,
1011 }
1012
1013
1014
1015 if(tbl_ref == NULL) {
1016 xsh_msg_error("Provide std sar catalog frame");
1017 return NULL;
1018 //cpl_table_dump(tbl_ref,0,3,stdout);
1019 }
1020
1021 for(i=0;i<next;i+=3) {
1022
1023
1024 //xsh_msg("Processing extracted spectrum extension %d",i);
1025 xsh_free_vector(&vec_ord);
1026 pobj=NULL;
1027 xsh_free_propertylist(&x_plist);
1028 check(vec_ord=cpl_vector_load(name_sci,i));
1029 check(x_plist=cpl_propertylist_load(name_sci,i));
1030 vec_size=cpl_vector_get_size(vec_ord);
1031
1032 //xsh_msg("naxis1=%d vec_ord nx=%d",naxis1,vec_size);
1033 check(pobj=cpl_vector_get_data(vec_ord));
1034 check(crval1=xsh_pfits_get_crval1(x_plist));
1035 check(cdelt1=xsh_pfits_get_cdelt1(x_plist));
1036 xsh_free_table(&tbl_ord);
1037 tbl_ord=cpl_table_new(vec_size);
1038 check(cpl_table_copy_structure(tbl_ord,obj_tab));
1039 check(cpl_table_fill_column_window_int(tbl_ord,COL_NAME_ORD_OBJ,0,vec_size,-1));
1040 check(cpl_table_fill_column_window_double(tbl_ord,COL_NAME_WAVE_OBJ,0,vec_size,-1));
1041 check(cpl_table_fill_column_window_double(tbl_ord,COL_NAME_INT_OBJ,0,vec_size,-1));
1042
1043 check(po=cpl_table_get_data_int(tbl_ord,COL_NAME_ORD_OBJ));
1044 check(pw=cpl_table_get_data_double(tbl_ord,COL_NAME_WAVE_OBJ));
1045 check(pf=cpl_table_get_data_double(tbl_ord,COL_NAME_INT_OBJ));
1046
1047
1048 //xsh_msg("crval1=%g cdelt1=%g",crval1,cdelt1);
1049 ord=i/3;
1050 for(j=0;j<vec_size;j++) {
1051 po[j]=ord;
1052 wav=crval1+cdelt1*j;
1053 pw[j]=wav;
1054 pf[j]=pobj[j];
1055 //xsh_msg("j=%d Flux[%g]=%g",j,pw[j],pobj[j]);
1056 }
1057
1058 //cpl_table_dump(tbl_ord,0,3,stdout);
1059 //sprintf(name_eff,"obj_ord%d.fits",ord);
1060 //check(cpl_table_save(tbl_ord,plist, x_plist,name_eff, CPL_IO_DEFAULT));
1061
1062 //sprintf(name_eff,"ref_ord%d.fits",ord);
1063 //check(cpl_table_save(tbl_ref,plist, x_plist,name_eff, CPL_IO_DEFAULT));
1064
1065 //check(cpl_table_save(tbl_atmext,plist, x_plist,"table_atm.fits", CPL_IO_DEFAULT));
1066
1067 xsh_free_table(&tbl_eff);
1068 check(tbl_eff=xsh_utils_efficiency_internal(tbl_ord,tbl_atmext,tbl_ref,
1069 exptime,airmass,aimprim,gain,
1070 biny,nm2nm,
1078 &ntot,&nclip));
1079
1080
1081 xsh_free_table(&tbl_ord);
1082 //cpl_table_dump(tbl_eff,1,2,stdout);
1083 nrow=cpl_table_get_nrow(tbl_eff);
1084 if(nrow>0) {
1085 check(emax=cpl_table_get_column_max(tbl_eff,"EFF"));
1086 } else {
1087 emax=-999;
1088 }
1089 sprintf(key_name,"%s%2d", XSH_QC_EFF_PEAK_ORD,ord);
1090 cpl_propertylist_append_double(plist,key_name,emax);
1091 cpl_propertylist_set_comment(plist,key_name,"Peak efficiency");
1092 if(nrow>0) {
1093 check(emed=cpl_table_get_column_median(tbl_eff,"EFF"));
1094 } else {
1095 emed=-999;
1096 }
1097 sprintf(key_name,"%s%2d", XSH_QC_EFF_MED_ORD,ord);
1098 cpl_propertylist_append_double(plist,key_name,emed);
1099 cpl_propertylist_set_comment(plist,key_name,"Median efficiency");
1100
1101
1102 neff_tot+=ntot;
1103 //xsh_msg("neff_tot=%d vec_size=%d",neff_tot,vec_size);
1104 if(ord==0) {
1105 tot_eff=cpl_table_duplicate(tbl_eff);
1106 } else {
1107 cpl_table_insert(tot_eff,tbl_eff,neff_tot);
1108 }
1109 nclip_tot+=nclip;
1110
1111 //sprintf(name_eff,"eff_ord%d.fits",ord);
1112 //check(cpl_table_save(tbl_eff,plist, x_plist,name_eff, CPL_IO_DEFAULT));
1113 }
1114
1115 //abort();
1116
1117 xsh_msg("nclip_tot=%d",nclip_tot);
1118 HIGH_ABS_REGION * phigh=NULL;
1119 phigh=xsh_fill_high_abs_regions(instrument,high_abs_win);
1120
1122
1123 sprintf(tag,"EFFICIENCY_%s_%s",xsh_instrument_mode_tostring(instrument),
1125 sprintf(fname,"%s.fits",tag);
1126
1127 nrow=cpl_table_get_nrow(tot_eff);
1128 fclip=nclip_tot;
1129
1130 if(neff_tot > 0) {
1131 fclip/=neff_tot;
1132 } else {
1133 fclip=-999;
1134 }
1135
1136 xsh_msg("nclip_tot =%d neff_tot=%d fclip=%g",nclip_tot,neff_tot,fclip);
1137 xsh_pfits_set_qc_eff_fclip(plist,fclip);
1138 xsh_pfits_set_qc_eff_nclip(plist,nclip_tot);
1139
1140 //cpl_propertylist_dump(plist,stdout);
1141 check(cpl_table_save(tot_eff,plist, x_plist,fname, CPL_IO_DEFAULT));
1142
1143 xsh_free_table(&obj_tab);
1144
1145 check(frm_eff=xsh_frame_product(fname,tag,CPL_FRAME_TYPE_TABLE,
1146 CPL_FRAME_GROUP_CALIB,
1147 CPL_FRAME_LEVEL_FINAL));
1148
1149
1150
1151 cleanup:
1152
1153 //if(phigh!=NULL) cpl_free(phigh);
1154
1155
1156 xsh_free_table(&tot_eff);
1157 xsh_free_table(&tbl_ref);
1158 xsh_free_table(&tbl_atmext);
1159 xsh_free_table(&obj_tab);
1160 xsh_free_table(&tbl_ord);
1161 xsh_free_table(&tbl_eff);
1162
1163 xsh_free_image(&ima_sci);
1165 xsh_free_vector(&vec_ord);
1166 xsh_free_image(&ima_obj);
1167 xsh_free_propertylist(&plist);
1168 xsh_free_propertylist(&x_plist);
1169
1170 return frm_eff;
1171
1172}
1173
1174
1175#if 0
1176static double slaAirmas ( double zd )
1177/*
1178** - - - - - - - - - -
1179** s l a A i r m a s
1180** - - - - - - - - - -
1181**
1182** Air mass at given zenith distance.
1183**
1184** (double precision)
1185**
1186** Given:
1187** zd d observed zenith distance (radians)
1188**
1189** The result is an estimate of the air mass, in units of that
1190** at the zenith.
1191**
1192** Notes:
1193**
1194** 1) The "observed" zenith distance referred to above means "as
1195** affected by refraction".
1196**
1197** 2) Uses Hardie's (1962) polynomial fit to Bemporad's data for
1198** the relative air mass, X, in units of thickness at the zenith
1199** as tabulated by Schoenberg (1929). This is adequate for all
1200** normal needs as it is accurate to better than 0.1% up to X =
1201** 6.8 and better than 1% up to X = 10. Bemporad's tabulated
1202** values are unlikely to be trustworthy to such accuracy
1203** because of variations in density, pressure and other
1204** conditions in the atmosphere from those assumed in his work.
1205**
1206** 3) The sign of the ZD is ignored.
1207**
1208** 4) At zenith distances greater than about ZD = 87 degrees the
1209** air mass is held constant to avoid arithmetic overflows.
1210**
1211** References:
1212** Hardie, R.H., 1962, in "Astronomical Techniques"
1213** ed. W.A. Hiltner, University of Chicago Press, p180.
1214** Schoenberg, E., 1929, Hdb. d. Ap.,
1215** Berlin, Julius Springer, 2, 268.
1216**
1217** Adapted from original Fortran code by P.W.Hill, St Andrews.
1218**
1219** Last revision: 5 October 1994
1220**
1221** Copyright P.T.Wallace. All rights reserved.
1222*/
1223{
1224 double w, seczm1;
1225
1226 w = fabs ( zd );
1227 /* replace
1228 seczm1 = 1.0 / ( cos ( gmin ( 1.52,w ) ) ) - 1.0;
1229 by :
1230 w = (1.52 > w) ? w : 1.52;
1231 seczm1 = 1.0 / ( cos ( w ) ) - 1.0;
1232 */
1233
1234 w = (1.52 > w) ? w : 1.52;
1235 seczm1 = 1.0 / ( cos ( w ) ) - 1.0;
1236 return 1.0 + seczm1 * ( 0.9981833
1237 - seczm1 * ( 0.002875 + 0.0008083 * seczm1 ) );
1238}
1239#endif
1240/*
1241 LST
1242 =========================
1243 * LST (sidereal time) keyword is retrieved from the TCS when the
1244 EXPSTRT is executed on the instrument side. Such commands is executed
1245 whenever and exposure is started on the instrument.
1246 * The keyword is added by the tif module, in methods
1247 tifFitsLib::tifGetFitsStart()
1248 * The value of the keyword, is the current (instantaneous) value found
1249 in telescope's DB:
1250 "trk:data:times.lst". This value is computed given the
1251 UTC from the time server:
1252*/
1253
1254
1255
1256void slaAoppat ( double date, double aoprms[14] )
1257/*
1258** - - - - - - - - - -
1259** s l a A o p p a t
1260** - - - - - - - - - -
1261**
1262** Recompute the sidereal time in the apparent to observed place
1263** star-independent parameter block.
1264**
1265** Given:
1266** date double UTC date/time (Modified Julian Date, JD-2400000.5)
1267** (see slaAoppa source for comments on leap seconds)
1268** aoprms double[14] star-independent apparent-to-observed parameters
1269**
1270** (0-11) not required
1271** (12) longitude + eqn of equinoxes + sidereal dut
1272** (13) not required
1273**
1274** Returned:
1275** aoprms double[14] star-independent apparent-to-observed parameters:
1276**
1277** (0-12) not changed
1278** (13) local apparent sidereal time (radians)
1279**
1280** For more information, see slaAoppa.
1281**
1282** Called: slaGmst
1283**
1284** Last revision: 31 October 1993
1285**
1286** Copyright P.T.Wallace. All rights reserved.
1287*/
1288{
1289 /* Temporarily commented out
1290 aoprms[13] = slaGmst ( date ) + aoprms[12];
1291 */
1292}
1293
1294static double
1295xsh_utils_get_airm(cpl_propertylist* plist)
1296{
1297
1298 double airmass=0;
1299 /* an unknown */
1300 double azimuth=0;
1301 /* an unknown */
1302 double altitude=0;
1303 /* From FITS header */
1304 double declination=0;
1305 double RA=0;
1306
1307 double hour_angle=0;
1308
1309 /* to get wich UT is used, that affects the latitude value */
1310 const char* telescope_id=NULL;
1311 int tel_id=0;
1312 double latitude=0;
1313 double lat_sig=-1;
1314 double lat_deg=24;
1315 double lat_min=37;
1316 double lat_sec=30; /* approximative value refined later according to UT# */
1317 int len=0;
1318 //double exptime_sec=0;
1319 //double exptime_hms=0;
1320 //double exptime_deg=0;
1321 double LST=0; /*Local Sidereal Time */
1322
1323 declination=xsh_pfits_get_dec(plist);
1324 RA=xsh_pfits_get_ra(plist);
1325 LST=xsh_pfits_get_lst(plist);
1326 hour_angle=LST-RA;
1327
1328 //xsh_msg("declination=%g",declination);
1329 //xsh_msg("right ascension=%g",RA);
1330 //xsh_msg("LST=%g",LST);
1331
1332 //xsh_msg("hour_angle=%g [hms]",hour_angle);
1333 //xsh_msg("hour_angle=%g [deg]",xsh_hms2deg(hour_angle));
1334
1335 telescope_id=xsh_pfits_get_telescop(plist);
1336 //exptime_sec=xsh_pfits_get_exptime(plist);
1337 //exptime_hms=exptime_sec/3600;
1338 //xsh_msg("exptime=%g [sec]",exptime_sec);
1339 //xsh_msg("exptime=%g [hms]",exptime_hms);
1340 //exptime_deg=xsh_hms2deg(exptime_hms);
1341 //xsh_msg("exptime=%g [deg]",exptime_deg);
1342
1343 //hour_angle+=exptime_deg;
1344
1345 //xsh_msg("telescope_id=%s",telescope_id);
1346 len=strlen(telescope_id);
1347 /* 48 is bitcode of 0 */
1348 tel_id=telescope_id[len-1]-48;
1349
1350 //xsh_msg("tel_id=%d",tel_id);
1351
1352 switch ( tel_id ) {
1353 case 1: lat_sec=33.117;
1354 break;
1355 case 2: lat_sec=31.465;
1356 break;
1357 case 3: lat_sec=30.300;
1358 break;
1359 case 4: lat_sec=31.000;
1360 break;
1361 }
1362 //xsh_msg("lat_sec=%g",lat_sec);
1363
1364 latitude=lat_sig*(lat_deg+lat_min/60+lat_sec/3600);
1365 //xsh_msg("latitude=%18.15g",latitude);
1366
1367 double c1=0;
1368 double c2=0;
1369 double c3=0;
1370
1371 c1=sin(latitude)*sin(declination)+cos(latitude)*cos(declination)*cos(hour_angle);
1372 c2=cos(latitude)*sin(declination)-sin(latitude)*cos(declination)*cos(hour_angle);
1373 c3=-cos(latitude)*sin(hour_angle);
1374
1375 azimuth=atan(c3/c2);
1376 altitude=atan(c1/c2*cos(azimuth));
1377 //xsh_msg("azimuth=%g",azimuth);
1378 xsh_msg("altitude=%g",altitude);
1379 //double zenith_distance=90-altitude;
1380 //xsh_msg("zenit_distance=%g",zenith_distance);
1381 /* Pickering 2002: X=1/sin(h+2444/(165+47h^1.1))
1382 see Wikipedia Air mass, interpolation formulae
1383 */
1384 airmass=1/sin(altitude+244/(165+47*pow(altitude,1.1)));
1385 //xsh_msg("airmass med=%g",airmass);
1386
1387 return airmass;
1388
1389
1390}
1391
1392
1393double
1395{
1396
1397 int i=0;
1398 int size=0;
1399 cpl_frame* frm=NULL;
1400 cpl_propertylist* plist=NULL;
1401 const char* name=NULL;
1402 double airm_start=0;
1403 double airm_end=0;
1404 double airm_med=0;
1405 double airm_eff=0;
1406
1407 double exptime=0;
1408 double area=0;
1409 double sum_area=0;
1410 double sum_time=0;
1411
1412 XSH_ASSURE_NOT_NULL_MSG(raws,"raws frameset null pointer");
1413
1414 size=cpl_frameset_get_size(raws);
1415 for(i=0;i<size;i++){
1416
1417 frm=cpl_frameset_get_frame(raws,i);
1418 name=cpl_frame_get_filename(frm);
1419 plist=cpl_propertylist_load(name,0);
1420
1421 airm_med=xsh_utils_get_airm(plist);
1422 airm_start=xsh_pfits_get_airm_start(plist);
1423 airm_end=xsh_pfits_get_airm_end(plist);
1425
1426 area=(exptime/6)*(airm_start+4*airm_med+airm_end);
1427
1428 //xsh_msg("area=%g",area);
1429 //xsh_msg("exptime=%g",exptime);
1430 sum_area+=area;
1431 sum_time+=exptime;
1432 xsh_free_propertylist(&plist);
1433
1434 }
1435 //xsh_msg("sum exptime=%g",sum_time);
1436 airm_eff=sum_area/sum_time;
1437
1438 cleanup:
1439 xsh_free_propertylist(&plist);
1440 xsh_msg("airmass eff=%g",airm_eff);
1441
1442 return airm_eff;
1443
1444}
1445
1446double
1447xsh_utils_compute_airm(cpl_frameset* raws)
1448{
1449
1450 cpl_frame* frm=NULL;
1451 cpl_propertylist* plist=NULL;
1452 const char* name=NULL;
1453 int i=0;
1454 int size=0;
1455 double airm=0;
1456 double airm_eff=0;
1457 double exptime=0;
1458 double time_tot=0;
1459 double* pairm=NULL;
1460 double* pprod=NULL;
1461 double* ptime=NULL;
1462
1463 cpl_array* time_array=NULL;
1464 cpl_array* airm_array=NULL;
1465 //cpl_array* prod_array=NULL;
1466
1467 XSH_ASSURE_NOT_NULL_MSG(raws,"raws frameset null pointer");
1468
1469 size=cpl_frameset_get_size(raws);
1470 airm_array=cpl_array_new(size,CPL_TYPE_DOUBLE);
1471 time_array=cpl_array_new(size,CPL_TYPE_DOUBLE);
1472 //prod_array=cpl_array_new(size,CPL_TYPE_DOUBLE);
1473 pairm=cpl_array_get_data_double(airm_array);
1474 ptime=cpl_array_get_data_double(time_array);
1475 pprod=cpl_array_get_data_double(time_array);
1476
1477 if(size>2) {
1478 for(i=0;i<size;i++){
1479
1480 frm=cpl_frameset_get_frame(raws,i);
1481 name=cpl_frame_get_filename(frm);
1482 plist=cpl_propertylist_load(name,0);
1483 airm=xsh_pfits_get_airm_mean(plist);
1485 ptime[i]=exptime;
1486 pairm[i]=airm;
1487 pprod[i]=airm*exptime;
1488
1489 }
1490 airm_eff=(pprod[0]+pprod[size-1])/(ptime[0]+ptime[size-1]);
1491
1492 } else if (size == 2) {
1493 frm=cpl_frameset_get_frame(raws,0);
1494 name=cpl_frame_get_filename(frm);
1495 plist=cpl_propertylist_load(name,0);
1496
1497 airm=xsh_pfits_get_airm_mean(plist);
1499
1500 airm_eff=airm*exptime;
1501 time_tot=exptime;
1502
1503 frm=cpl_frameset_get_frame(raws,1);
1504 name=cpl_frame_get_filename(frm);
1505 plist=cpl_propertylist_load(name,0);
1506
1507
1508 airm_eff+=airm*exptime;
1509 time_tot+=exptime;
1510
1511 airm_eff/=time_tot;
1512 } else {
1513
1514 frm=cpl_frameset_get_frame(raws,0);
1515 name=cpl_frame_get_filename(frm);
1516 plist=cpl_propertylist_load(name,0);
1517
1518 airm=xsh_pfits_get_airm_mean(plist);
1520
1521 airm_eff=airm;
1522 }
1523
1524 cleanup:
1525 return airm_eff;
1526
1527}
1528
1529
static double exptime
static xsh_instrument * instrument
int biny
#define XSH_ASSURE_NOT_NULL_MSG(pointer, msg)
Definition: xsh_error.h:103
#define check(COMMAND)
Definition: xsh_error.h:71
const char * xsh_instrument_mode_tostring(xsh_instrument *i)
Get the string associated with a mode.
const char * xsh_instrument_arm_tostring(xsh_instrument *i)
Get the string associated with an arm.
XSH_ARM xsh_instrument_get_arm(xsh_instrument *i)
Get an arm on instrument structure.
int size
#define xsh_msg_warning(...)
Print an warning message.
Definition: xsh_msg.h:88
#define xsh_msg_dbg_medium(...)
Definition: xsh_msg.h:44
#define xsh_msg_error(...)
Print an error message.
Definition: xsh_msg.h:62
#define xsh_msg(...)
Print a message on info level.
Definition: xsh_msg.h:121
double xsh_pfits_get_lst(const cpl_propertylist *plist)
find out the LST value
Definition: xsh_pfits.c:478
double xsh_pfits_get_conad(const cpl_propertylist *plist)
find out the CONAD value
Definition: xsh_pfits.c:692
double xsh_pfits_get_airm_end(const cpl_propertylist *plist)
find out the TEL AIRM END value
Definition: xsh_pfits.c:527
double xsh_pfits_get_dit(const cpl_propertylist *plist)
find out the DIT value
Definition: xsh_pfits.c:1405
double xsh_pfits_get_win1_dit1(const cpl_propertylist *plist)
find out the DET WIN1 DIT1 value
Definition: xsh_pfits.c:1385
double xsh_pfits_get_dec(const cpl_propertylist *plist)
Get the Right Ascension.
Definition: xsh_pfits.c:3500
const char * xsh_pfits_get_telescop(const cpl_propertylist *plist)
find out the TELESCOP value (telescope ID)
Definition: xsh_pfits.c:657
double xsh_pfits_get_cdelt1(const cpl_propertylist *plist)
find out the cdelt1
Definition: xsh_pfits.c:2196
int xsh_pfits_get_biny(const cpl_propertylist *plist)
find out the BINY value
Definition: xsh_pfits.c:306
double xsh_pfits_get_exptime(const cpl_propertylist *plist)
find out the exposure time
Definition: xsh_pfits.c:2254
double xsh_pfits_get_ra(const cpl_propertylist *plist)
Get the Right Ascension.
Definition: xsh_pfits.c:3479
double xsh_pfits_get_airm_mean(const cpl_propertylist *plist)
find out the mean airmass value
Definition: xsh_pfits.c:511
int xsh_pfits_get_naxis1(const cpl_propertylist *plist)
find out the NAXIS1 value
Definition: xsh_pfits.c:227
double xsh_pfits_get_crval1(const cpl_propertylist *plist)
find out the crval1
Definition: xsh_pfits.c:1907
double xsh_pfits_get_airm_start(const cpl_propertylist *plist)
find out the TEL AIRM START value
Definition: xsh_pfits.c:496
void xsh_pfits_set_qc_eff_fclip(cpl_propertylist *plist, double value)
void xsh_pfits_set_qc_eff_nclip(cpl_propertylist *plist, int value)
static cpl_vector * rec_profile
Definition: xsh_rectify.c:80
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
double xsh_table_interpolate(cpl_table *tbl, double wav, const char *colx, const char *coly)
Interpolate table columns.
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
#define XSH_ATMOS_EXT_LIST_COLNAME_OLD
#define XSH_ATMOS_EXT_LIST_COLNAME_K
@ XSH_ARM_AGC
@ XSH_ARM_UNDEFINED
@ XSH_ARM_UVB
@ XSH_ARM_NIR
@ XSH_ARM_VIS
double kappa
Definition: xsh_detmon_lg.c:81
cpl_frame * xsh_frame_product(const char *fname, const char *tag, cpl_frame_type type, cpl_frame_group group, cpl_frame_level level)
Creates a frame with given characteristics.
Definition: xsh_dfs.c:930
cpl_frame * xsh_find_frame_with_tag(cpl_frameset *frames, const char *tag, xsh_instrument *instr)
Find frame with a given tag.
Definition: xsh_dfs.c:3347
#define XSH_FLUX_STD_CAT
Definition: xsh_dfs.h:272
#define XSH_STD_FLUX_SLIT_STARE_ORDER1D
Definition: xsh_dfs.h:309
#define XSH_EXTCOEFF_TAB
Definition: xsh_dfs.h:277
#define XSH_FLUX_STD_TAB
Definition: xsh_dfs.h:266
HIGH_ABS_REGION * xsh_fill_high_abs_regions(xsh_instrument *instrument, cpl_frame *high_abs_frame)
static const double STAR_MATCH_DEPSILON
#define XSH_QC_EFF_MED_ORD
Definition: xsh_pfits_qc.h:49
#define XSH_QC_EFF_PEAK_ORD
Definition: xsh_pfits_qc.h:50
void star_index_delete(star_index *pindex)
star_index * star_index_load(const char *fits_file)
cpl_table * star_index_get(star_index *pindex, double RA, double DEC, double RA_EPS, double DEC_EPS, const char **pstar_name)
#define TEL_AREA
Definition: xsh_util_afc.c:63
double xsh_utils_compute_airm(cpl_frameset *raws)
static const char COL_NAME_COR[]
void xsh_load_ref_table(cpl_frameset *frames, double dRA, double dDEC, double EPSILON, xsh_instrument *instrument, cpl_table **pptable)
load reference table
static const char COL_NAME_EXT[]
double xsh_utils_compute_airm_eff(cpl_frameset *raws)
static void xsh_frame_sci_get_gain_airmass_exptime_naxis1_biny(cpl_frame *frm_sci, xsh_instrument *instrument, double *gain, double *airmass, double *exptime, int *naxis1, int *biny)
static char FRM_EXTCOEFF_TAB[]
static const char COL_NAME_ABS_ATMDISP[]
static cpl_error_code xsh_get_std_obs_values(cpl_propertylist *plist, double *exptime, double *airmass, double *dRA, double *dDEC)
get STD star observation exptime, airmass, RA, DEC
cpl_error_code xsh_efficiency_add_high_abs_regions(cpl_table **eff, HIGH_ABS_REGION *phigh)
static const char COL_NAME_WAVELENGTH[]
cpl_error_code xsh_rv_ref_wave_init(xsh_std_star_id std_star_id, XSH_ARM arm, xsh_rv_ref_wave_param *w)
static const char COL_NAME_BINWIDTH_REF[]
cpl_frame * xsh_utils_efficiency(cpl_frameset *frames, double dGain, double dEpsilon, double aimprim, xsh_instrument *inst, const char *col_name_atm_wave, const char *col_name_atm_abs, const char *col_name_ref_wave, const char *col_name_ref_flux, const char *col_name_ref_bin, const char *col_name_obj_wave, const char *col_name_obj_flux)
Compute efficiency.
static const char COL_NAME_ORD_OBJ[]
static const char COL_NAME_WAVE_ATMDISP[]
static const char COL_NAME_EPHOT[]
static double xsh_utils_get_airm(cpl_propertylist *plist)
static const char COL_NAME_FLUX_REF[]
static const char COL_NAME_WAVE_REF[]
void slaAoppat(double date, double aoprms[14])
void xsh_rv_ref_wave_param_destroy(xsh_rv_ref_wave_param *p)
cpl_table * xsh_utils_efficiency_internal(cpl_table *tbl_obj_spectrum, cpl_table *tbl_atmext, cpl_table *tbl_ref, double exptime, double airmass, double aimprim, double gain, int biny, double src2ref_wave_sampling, const char *col_name_atm_wave, const char *col_name_atm_abs, const char *col_name_ref_wave, const char *col_name_ref_flux, const char *col_name_ref_bin, const char *col_name_obj_wave, const char *col_name_obj_flux, int *ntot, int *nclip)
Compute efficiency.
void xsh_frame_sci_get_ra_dec_airmass(cpl_frame *frm_sci, double *ra, double *dec, double *airmass)
get RA, DEC, airmass (mean) of a frame
static double * xsh_create_column_double(cpl_table *tbl, const char *col_name, int nrow)
static const char COL_NAME_WAVE_OBJ[]
static const char COL_NAME_INT_OBJ[]
cpl_frame * xsh_efficiency_compute(cpl_frame *frm_sci, cpl_frame *frm_cat, cpl_frame *frm_atmext, cpl_frame *high_abs_win, xsh_instrument *instrument)
computes efficiency
cpl_frame * xsh_catalog_extract_spectrum_frame(cpl_frame *frm_cat, cpl_frame *frm_sci)
extract spectrum
static const char COL_NAME_SRC_COR[]
xsh_rv_ref_wave_param * xsh_rv_ref_wave_param_create(void)
static int xsh_column_to_double(cpl_table *ptable, const char *column)
cpl_error_code xsh_parse_catalog_std_stars(cpl_frame *cat, double dRA, double dDEC, double EPSILON, cpl_table **pptable, xsh_std_star_id *std_star_id)
parse referece std stars catalog
static const char COL_NAME_REF[]
static const char COL_NAME_HIGH_ABS[]
static const char COL_NAME_SRC_EFF[]
xsh_std_star_id
@ XSH_LTT3218
@ XSH_GD71
@ XSH_LTT7987
@ XSH_Feige110
@ XSH_BD17
@ XSH_EG274
@ XSH_GD153