visir_spc_photom.c

00001 /* $Id: visir_spc_photom.c,v 1.53 2009/03/06 13:42:03 llundin Exp $
00002  *
00003  * This file is part of the VISIR Pipeline
00004  * Copyright (C) 2002,2003 European Southern Observatory
00005  *
00006  * This program is free software; you can redistribute it and/or modify
00007  * it under the terms of the GNU General Public License as published by
00008  * the Free Software Foundation; either version 2 of the License, or
00009  * (at your option) any later version.
00010  *
00011  * This program is distributed in the hope that it will be useful,
00012  * but WITHOUT ANY WARRANTY; without even the implied warranty of
00013  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00014  * GNU General Public License for more details.
00015  *
00016  * You should have received a copy of the GNU General Public License
00017  * along with this program; if not, write to the Free Software
00018  * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02111-1307  USA
00019  */
00020 
00021 /*
00022  * $Author: llundin $
00023  * $Date: 2009/03/06 13:42:03 $
00024  * $Revision: 1.53 $
00025  * $Name: visir-3_2_2 $
00026  */
00027 
00028 #ifdef HAVE_CONFIG_H
00029 #include <config.h>
00030 #endif
00031 
00032 /*-----------------------------------------------------------------------------
00033                                        Includes
00034  -----------------------------------------------------------------------------*/
00035 
00036 #include <string.h>
00037 #include <math.h>
00038 #include <assert.h>
00039 
00040 #include "irplib_framelist.h"
00041 
00042 #include "visir_spc_photom.h"
00043 
00044 #include "visir_parameter.h"
00045 #include "visir_utils.h"
00046 #include "visir_pfits.h"
00047 #include "visir_spc_distortion.h"
00048 
00049 /*-----------------------------------------------------------------------------
00050                             Static variables
00051  -----------------------------------------------------------------------------*/
00052 
00053 static const int fits_strlen = IRPLIB_FITS_STRLEN;
00054 
00055 /*-----------------------------------------------------------------------------
00056                           Private Function prototypes
00057  -----------------------------------------------------------------------------*/
00058 
00059 static cpl_bivector * visir_spc_phot_model_from_cat(const char *, double,
00060                                                     double);
00061 static int visir_spc_phot_model_rebin(const cpl_bivector *, cpl_table *);
00062 
00063 static cpl_error_code visir_spc_sens_stat(const visir_spc_phot_config *,
00064                                           double *, double *,
00065                                           const cpl_table *);
00066 
00067 static cpl_error_code visir_spc_phot_qc(cpl_propertylist *,
00068                                         const char *,
00069                                         double, double, double);
00070 
00071 /*----------------------------------------------------------------------------*/
00075 /*----------------------------------------------------------------------------*/
00076 
00079 /*----------------------------------------------------------------------------*/
00100 /*----------------------------------------------------------------------------*/
00101 cpl_table * visir_spc_phot_sensit(
00102         cpl_image    ** combined,
00103         const irplib_framelist * rawframes,
00104         const char   * recipename,
00105         const cpl_parameterlist * parlist,
00106         const visir_spc_phot_config * pconfig,
00107         const char   * star_cat,
00108         const char   * spc_cal_lines,
00109         const char   * spc_cal_qeff,
00110         cpl_image   ** pweight2d,
00111         cpl_propertylist * qclist,
00112         cpl_boolean    do_ech,
00113         double         wlen,
00114         double         slitw,
00115         double         temp,
00116         double         fwhm,
00117         const visir_spc_resol resol)
00118 {
00119     cpl_errorstate cleanstate = cpl_errorstate_get();
00120     const cpl_propertylist * plist;
00121     double                  exptime;
00122     double                  dit, ndit;
00123     int                     ncycles, nnod;
00124     int                     icol1, icol2;
00125     int                     jcol1, jcol2;
00126     double                  ra, dec;
00127     cpl_imagelist       *   hcycle = NULL;
00128     cpl_image           *   imhcycle = NULL;
00129     cpl_image           *   comorder = NULL;
00130     cpl_image           *   order = NULL;
00131     const char          *   star_name;
00132     cpl_bivector        *   flux_model = NULL;
00133     cpl_table           *   spc_table = NULL;
00134     double                  factor;
00135     double              *   ext_array;
00136     double              *   spc_array;
00137     double              *   err_array;
00138     double                  sens_mean = -1.0; /* Avoid (false) uninit warning */
00139     double                  sens_stdev = -1.0; /* Ditto */
00140     int                     i;
00141     
00142 
00143     *pweight2d = NULL;
00144 
00145     skip_if (0);
00146 
00147     skip_if (rawframes == NULL);
00148 
00149     jcol1 = visir_parameterlist_get_int(parlist, recipename,
00150                                         VISIR_PARAM_REJLEFT);
00151     jcol2 = visir_parameterlist_get_int(parlist, recipename,
00152                                         VISIR_PARAM_REJRIGHT);
00153 
00154     skip_if (0);
00155 
00156     /* DIT must be present every where */  
00157     skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_DOUBLE_DIT,
00158                                       CPL_TYPE_DOUBLE, CPL_FALSE, 0.0));
00159     if(irplib_framelist_contains(rawframes, VISIR_PFITS_DOUBLE_DIT,
00160                                  CPL_TYPE_DOUBLE, CPL_TRUE, 1e-5)) {
00161         /* Allow 0.1 ms difference - or warn */
00162         /* FIXME: The recipe does not properly handle non-uniform DITSs */  
00163         visir_error_reset("DIT differs by more than %g", 1e-5);
00164     }
00165 
00166     /* FIXME: Verify the angular distance */
00167     skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_DOUBLE_RA,
00168                                       CPL_TYPE_DOUBLE, CPL_FALSE, 0.0));
00169 
00170     /* FIXME: Allow 1 degree difference */
00171     skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_DOUBLE_DEC,
00172                                       CPL_TYPE_DOUBLE, CPL_TRUE, 1.0));
00173 
00174     skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_INT_CHOP_NCYCLES,
00175                                       CPL_TYPE_INT, CPL_TRUE, 0.0));
00176 
00177     skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_INT_NDIT,
00178                                       CPL_TYPE_INT, CPL_TRUE, 0.0));
00179 
00180     if (irplib_framelist_contains(rawframes, VISIR_PFITS_STRING_STARNAME,
00181                                   CPL_TYPE_STRING, CPL_TRUE, 0.0)) {
00182         visir_error_reset("Rawframe(s) missing standard star name");
00183     }
00184 
00185 
00186     /* Get the first frame */
00187     plist = irplib_framelist_get_propertylist_const(rawframes, 0);
00188 
00189     skip_if (0);
00190 
00191     /* Get the total exposure time */
00192     /* DIT */
00193     dit = visir_pfits_get_dit(plist);
00194     /* NDIT */
00195     ndit = visir_pfits_get_ndit(plist);
00196     /* NNOD */
00197     nnod = irplib_framelist_get_size(rawframes);
00198     /* Number of chopping cycles */
00199     ncycles = visir_pfits_get_chop_ncycles(plist);
00200     /* exptime *= 2 because of chopping */
00201     exptime = 2 * dit * ndit * nnod * ncycles;
00202 
00203     skip_if (0);
00204 
00205     /* Get the standard star name */
00206     star_name = visir_pfits_get_starname(plist);
00207     if (star_name == NULL) visir_error_reset("Could not get standard star name");
00208     
00209     /* Get RA / DEC */
00210     ra  = visir_pfits_get_ra(plist);
00211     dec = visir_pfits_get_dec(plist);
00212     
00213     skip_if (0);
00214 
00215     /* Get the HCYCLE frame */
00216     hcycle = visir_load_imagelist(rawframes, 0, CPL_FALSE);
00217     skip_if (0);
00218 
00219     skip_if (pconfig->do_fixcombi && visir_spc_det_fix(combined, 1, CPL_TRUE,
00220                                wlen, resol,
00221                                pconfig->phi,
00222                                pconfig->ksi,
00223                                pconfig->eps,
00224                                pconfig->delta,
00225                                pconfig->plot));
00226 
00227 
00228     imhcycle = cpl_imagelist_get(hcycle, 0);    
00229     skip_if (0);
00230    
00231     skip_if (visir_spc_det_fix(&imhcycle, 1, CPL_FALSE,
00232                                wlen, resol,
00233                                pconfig->phi,
00234                                pconfig->ksi,
00235                                pconfig->eps,
00236                                pconfig->delta,
00237                                pconfig->plot));
00238 
00239     if (do_ech) {
00240         skip_if (visir_spc_echelle_limit(&icol1, &icol2, wlen,
00241                                          pconfig->orderoffset,
00242                                          1, cpl_image_get_size_y(*combined)));
00243         skip_if (visir_qc_append_background(qclist, rawframes, icol1, icol2));
00244 
00245     } else {
00246         icol1 = 1;
00247         icol2 = cpl_image_get_size_x(imhcycle);
00248 
00249         skip_if (visir_qc_append_background(qclist, rawframes, 0, 0));
00250     }
00251 
00252     if (do_ech) {
00253         if (jcol1 != 0) {
00254             cpl_msg_info(cpl_func, "Ignoring %d leftmost columns from %d to %d",
00255                          jcol1, icol1, icol1 + jcol1);
00256             icol1 += jcol1;
00257         }
00258         if (jcol2 != 0) {
00259             cpl_msg_info(cpl_func, "Ignoring %d rightmost columns from %d to %d",
00260                          jcol2, icol2 - jcol2, icol2);
00261             icol2 -= jcol2;
00262         }
00263     } else {
00264         if (jcol1 != 0) {
00265             cpl_msg_info(cpl_func, "Ignoring %d leftmost columns", jcol1);
00266             icol1 += jcol1;
00267         }
00268         if (jcol2 != 0) {
00269             cpl_msg_info(cpl_func, "Ignoring %d rightmost columns", jcol2);
00270             icol2 -= jcol2;
00271         }
00272     }
00273 
00274     if (icol1 != 1 || icol2 != cpl_image_get_size_x(imhcycle)) {
00275         order = visir_spc_column_extract(imhcycle, icol1, icol2,
00276                                           pconfig->plot);
00277         skip_if (0);
00278 
00279         comorder = visir_spc_column_extract(*combined, icol1, icol2,
00280                                              pconfig->plot);
00281         skip_if (0);
00282 
00283     } else {
00284         order = imhcycle;
00285         comorder = *combined;
00286     }
00287 
00288     /* Extract the spectrum */
00289     skip_if (visir_spc_extract_wcal(comorder,
00290                                     order, wlen, slitw, temp, fwhm,
00291                                     resol, pconfig->orderoffset,
00292                                     spc_cal_lines, spc_cal_qeff,
00293                                     &spc_table, pweight2d,
00294                                     qclist,
00295                                     pconfig->plot));
00296       
00297     cpl_imagelist_delete(hcycle);
00298     hcycle = NULL;
00299 
00300     /* Get the flux model from the catalog */
00301     if ((flux_model=visir_spc_phot_model_from_cat(star_cat, ra, dec)) == NULL) {
00302         cpl_msg_error(cpl_func, "Cannot retrieve the flux model from the cat.");
00303         skip_if(1);
00304     }
00305 
00306     /* Extract the interesting part of the flux model and rebin it */
00307     if (visir_spc_phot_model_rebin(flux_model, spc_table) == -1) {
00308         cpl_msg_error(cpl_func, "Cannot rebin the flux model");
00309         skip_if(1);
00310     }
00311     cpl_bivector_delete(flux_model);
00312     flux_model = NULL;
00313 
00314     /* Compute the sensitivity */
00315     /* Sensitivity = model * error * 10 * sqrt(DIT*NDIT*NFILES*NCHOP*2) 
00316         / extracted */
00317     cpl_table_duplicate_column(spc_table, "SENSITIVITY", spc_table, 
00318             "STD_STAR_MODEL");
00319     cpl_table_set_column_unit(spc_table, "SENSITIVITY", "mJy");
00320 
00321     ext_array = cpl_table_get_data_double(spc_table, "SENSITIVITY");
00322     spc_array = cpl_table_get_data_double(spc_table, "SPC_EXTRACTED");
00323     err_array = cpl_table_get_data_double(spc_table, "SPC_ERROR");
00324 
00325     factor = sqrt(exptime)/6.0; /* sqrt(exptime/3600) * 10; */
00326     for (i = 0; i < cpl_table_get_nrow(spc_table); i++) {
00327         ext_array[i] *= factor * err_array[i];
00328         if (ext_array[i] < 0 || spc_array[i] <= 0) {
00329             cpl_msg_warning(cpl_func, "Setting non-physical sensitivity in row %d "
00330                             "to 0: %g/%g", i+1, ext_array[i], spc_array[i]);
00331             ext_array[i] = 0;
00332         } else {
00333             ext_array[i] /= spc_array[i];
00334         }
00335     }
00336 
00337     skip_if(visir_spc_sens_stat(pconfig, &sens_mean, &sens_stdev, spc_table));
00338 
00339     /* This column is not part of the product */
00340     skip_if (cpl_table_erase_column(spc_table, "SPC_EMISSIVITY"));
00341 
00342     skip_if (visir_spc_phot_qc(qclist, star_name, exptime, sens_mean,
00343                                sens_stdev));
00344 
00345     /* Plot the results */
00346     if (pconfig->plot) {
00347         visir_table_plot("", "t 'Extracted spectrum' w lines", "", spc_table,
00348                          "WLEN", "SPC_EXTRACTED");
00349         visir_table_plot("", "t 'Extracted spectrum error' w lines", "", 
00350                          spc_table, "WLEN", "SPC_ERROR");
00351         visir_table_plot("", "t 'Standard star model' w lines", "", spc_table, 
00352                          "WLEN", "STD_STAR_MODEL");
00353         visir_table_plot("set grid;", "t 'Sensitivity (mJy)' w lines", "",
00354                          spc_table, "WLEN", "SENSITIVITY");
00355         visir_image_plot("", "t 'The weight map'", "", *pweight2d);
00356     }
00357 
00358     end_skip;
00359 
00360     cpl_imagelist_delete(hcycle);
00361     cpl_bivector_delete(flux_model);
00362 
00363     if (order != imhcycle) cpl_image_delete(order);
00364     if (comorder != *combined) cpl_image_delete(comorder);
00365 
00366     if (cpl_error_get_code()) {
00367         cpl_table_delete(spc_table);
00368         spc_table = NULL;
00369         cpl_image_delete(*pweight2d);
00370         *pweight2d = NULL;
00371     }
00372 
00373     return spc_table;
00374 }
00375 
00378 /*----------------------------------------------------------------------------*/
00390 /*----------------------------------------------------------------------------*/
00391 static cpl_bivector * visir_spc_phot_model_from_cat(const char * star_cat,
00392                                                     double       ra,
00393                                                     double       dec)
00394 {
00395     const double      max_radius = VISIR_STAR_MAX_RADIUS;
00396     const cpl_table * catalog  = NULL;  /* Flux model for all stars */ 
00397     int               nb_stars;
00398     const double    * dd;
00399     const cpl_array * da;
00400     cpl_vector      * v_ra = NULL;
00401     cpl_vector      * v_dec = NULL;
00402     int               min_dist_ind;
00403     cpl_bivector    * model   = NULL;
00404     cpl_vector      * model_x = NULL;
00405     cpl_vector      * model_y = NULL;
00406     const double      conv_mJy = 3.33E8; /* FIXME: 1e8/3 ? */
00407     double            star_dist;
00408     int               nb_vals;
00409 
00410 
00411     bug_if (0);  
00412 
00413     /* Open the star catalog */
00414     catalog = cpl_table_load(star_cat, 1, 0);
00415     if (catalog == NULL) {
00416         cpl_msg_error(cpl_func, "Could not open the star catalog: %s",
00417                      star_cat ? star_cat : "<NULL>");
00418         skip_if(1);
00419     }
00420 
00421     nb_stars = cpl_table_get_nrow(catalog);
00422 
00423     skip_if(nb_stars < 1);
00424 
00425     /* Get the RA and DEC columns */
00426     dd  = cpl_table_get_data_double_const(catalog, "RA");
00427     skip_if (dd == NULL);
00428     v_ra = cpl_vector_wrap(nb_stars, (double*)dd); /* v_ra is _not_ modified */
00429     bug_if( v_ra  == NULL);
00430 
00431     dd  = cpl_table_get_data_double_const(catalog, "DEC");
00432     skip_if (dd == NULL);
00433     v_dec = cpl_vector_wrap(nb_stars, (double*)dd); /* v_dec is _not_ modified */
00434     bug_if( v_dec  == NULL);
00435 
00436     /* Find the star closest to the given ra, dec */
00437     min_dist_ind = visir_star_find(v_ra, v_dec, ra, dec, max_radius, &star_dist);
00438 
00439     skip_if (min_dist_ind < 0);
00440 
00441     cpl_msg_info(cpl_func, "The standard star closest to (RA,DEC)=(%g,%g) is "
00442                  "no. %d, '%s' at (RA,DEC)=(%g,%g) with the distance [degree]: "
00443                  "%g", ra, dec, 1+min_dist_ind,
00444                  cpl_table_get_string(catalog, "STARS", min_dist_ind),
00445                  cpl_table_get_double(catalog, "RA",    min_dist_ind, NULL),
00446                  cpl_table_get_double(catalog, "DEC",   min_dist_ind, NULL),
00447                  star_dist);
00448 
00449     /* Get the wavelengths of the star model */
00450     da = cpl_table_get_array(catalog, "WAVELENGTHS", min_dist_ind);
00451     skip_if (da == NULL);
00452     dd = cpl_array_get_data_double_const(da);
00453     skip_if (dd == NULL);
00454 
00455     nb_vals = cpl_array_get_size(da);
00456 
00457     model_x = cpl_vector_new(nb_vals);
00458     memcpy(cpl_vector_get_data(model_x), dd, nb_vals * sizeof(double));
00459 
00460     /* Get the corresponding fluxes of the star model */
00461     da = cpl_table_get_array(catalog, "MODEL_FLUX", min_dist_ind);
00462     skip_if (da == NULL);
00463     dd = cpl_array_get_data_double_const(da);
00464     skip_if (dd == NULL);
00465 
00466     skip_if (nb_vals != cpl_array_get_size(da));
00467 
00468     model_y = cpl_vector_new(nb_vals);
00469     memcpy(cpl_vector_get_data(model_y), dd, nb_vals * sizeof(double));
00470 
00471     /* Convert from W/m2/m to mJy with * 3.33 * 10^8 * lamda * lamda  */
00472     bug_if (cpl_vector_multiply_scalar(model_y, conv_mJy));
00473     bug_if (cpl_vector_multiply(model_y, model_x));
00474     bug_if (cpl_vector_multiply(model_y, model_x));
00475 
00476     /* Convert from microns to meters */
00477     bug_if (cpl_vector_multiply_scalar(model_x, 1e-6));
00478 
00479     model = cpl_bivector_wrap_vectors(model_x, model_y);
00480 
00481     bug_if (model == NULL);
00482 
00483     end_skip;
00484 
00485     (void)cpl_vector_unwrap(v_ra);
00486     (void)cpl_vector_unwrap(v_dec);
00487     cpl_table_delete((cpl_table*)catalog);
00488 
00489     if (cpl_error_get_code()) {
00490         /* At this point model is certain to be NULL */
00491         cpl_vector_delete(model_x);
00492         cpl_vector_delete(model_y);
00493     }   
00494 
00495     return model;
00496 }
00497 
00498 /*----------------------------------------------------------------------------*/
00504 /*----------------------------------------------------------------------------*/
00505 static int visir_spc_phot_model_rebin(
00506         const cpl_bivector * model,
00507         cpl_table          * spec_tab)
00508 {
00509     cpl_vector      *   bounds = NULL;
00510     cpl_vector      *   spec = NULL;
00511     double              bin_pos;
00512     const int           nrow = cpl_table_get_nrow(spec_tab);
00513     int                 i;
00514 
00515 
00516     skip_if (0);
00517 
00518     /* Create bounds */
00519     bounds = cpl_vector_new(nrow + 1);
00520     for (i=1 ; i<cpl_vector_get_size(bounds) - 1 ; i++) {
00521         bin_pos = (cpl_table_get(spec_tab, "WLEN", i-1, NULL) + 
00522                 cpl_table_get(spec_tab, "WLEN", i, NULL)) / 2.0;
00523         cpl_vector_set(bounds, i, bin_pos);
00524     }
00525     bin_pos = cpl_table_get(spec_tab, "WLEN", 0, NULL) - 
00526         ((cpl_table_get(spec_tab, "WLEN", 1, NULL) -
00527          cpl_table_get(spec_tab, "WLEN", 0, NULL)) / 2.0);
00528     cpl_vector_set(bounds, 0, bin_pos);
00529     bin_pos = 
00530         cpl_table_get(spec_tab, "WLEN", nrow-1, NULL) + 
00531         ((cpl_table_get(spec_tab,"WLEN", 1, NULL) -
00532          cpl_table_get(spec_tab, "WLEN", 0, NULL)) / 2.0);
00533     cpl_vector_set(bounds, cpl_vector_get_size(bounds)-1, bin_pos);
00534    
00535     /* Create the interpolated spectrum */
00536     spec = cpl_vector_new(nrow);
00537 
00538     /* Interpolate the spectrum */
00539     if (visir_vector_resample(spec, bounds, model)) {
00540         cpl_msg_error(cpl_func, "Cannot rebin the spectrum");
00541         skip_if(1);
00542     }
00543 
00544     /* Add the result in spec_tab */
00545     cpl_table_new_column(spec_tab, "STD_STAR_MODEL", CPL_TYPE_DOUBLE);
00546     cpl_table_set_column_unit(spec_tab, "STD_STAR_MODEL", "mJy");
00547     for (i=0 ; i<nrow ; i++) {
00548         cpl_table_set_double(spec_tab, "STD_STAR_MODEL", i, 
00549                 cpl_vector_get(spec, i));
00550     }
00551 
00552     end_skip;
00553 
00554     cpl_vector_delete(spec);
00555     cpl_vector_delete(bounds);
00556 
00557     return cpl_error_get_code();
00558 }
00559 
00560 /*----------------------------------------------------------------------------*/
00569 /*----------------------------------------------------------------------------*/
00570 static cpl_error_code visir_spc_sens_stat(const visir_spc_phot_config * pconfig,
00571                                           double * psens_mean,
00572                                           double * psens_stdev,
00573                                           const cpl_table * spc_table)
00574 {
00575 
00576     cpl_vector * sens_ignore = NULL;
00577     cpl_vector * sens_mini = NULL;
00578     double       emis_min, emis_max;
00579     const int    npix = cpl_table_get_nrow(spc_table);
00580     const int    npix_ignore = 3;
00581     int i;
00582     int iok;
00583 
00584     skip_if(0);
00585 
00586     skip_if(pconfig  == NULL);
00587 
00588     /* Compute mean and stdev of the sensitivity signal
00589        - ignore first and last npix_ignore elements */
00590     /* FIXME: Use cpl_table_select() */
00591     /* sens_ignore will _not_ be modified */
00592     sens_ignore = cpl_vector_wrap( npix - 2*npix_ignore, (double*)
00593                                cpl_table_get_data_double_const(spc_table,
00594                                                  "SENSITIVITY") + npix_ignore);
00595     skip_if(0);
00596 
00597     if (pconfig->plot)
00598         visir_vector_plot("set grid;","t 'Truncated Sensitivity (mJy)' w lines",
00599                           "", sens_ignore);
00600 
00601     sens_mini = cpl_vector_duplicate(sens_ignore);
00602     skip_if(0);
00603 
00604     emis_min = cpl_table_get_column_min(spc_table, "SPC_EMISSIVITY");
00605     emis_max = cpl_table_get_column_max(spc_table, "SPC_EMISSIVITY");
00606     skip_if(0);
00607 
00608     emis_max = emis_min + pconfig->emis_tol * (emis_max - emis_min);
00609 
00610     iok = 0;
00611     for (i=0; i < npix - 2*npix_ignore; i++) {
00612       const double emis = cpl_table_get(spc_table, "SPC_EMISSIVITY",
00613                                         i + npix_ignore, NULL);
00614 
00615       skip_if(0);
00616 
00617       if (emis > emis_max) continue;
00618 
00619       if (i > iok)
00620           skip_if(cpl_vector_set(sens_mini, iok, cpl_vector_get(sens_mini, i)));
00621 
00622       iok++;
00623 
00624     }
00625 
00626     assert( iok > 0);
00627 
00628     skip_if(cpl_vector_set_size(sens_mini, iok));
00629 
00630     *psens_mean =  cpl_vector_get_mean(sens_mini);
00631 
00632     if (iok == 1) {
00633         cpl_msg_warning(cpl_func, "Sensitivity computed on only 1 wavelength "
00634                      "with emissivity %f", emis_max);
00635         *psens_stdev = 0;
00636     } else {
00637         cpl_msg_info(cpl_func, "Sensitivity computed on %d wavelengths with "
00638                      "emissivity at most %f", iok, emis_max);
00639 
00640         *psens_stdev = cpl_vector_get_stdev(sens_mini);
00641 
00642     }
00643 
00644 
00645     end_skip;
00646 
00647     cpl_vector_unwrap(sens_ignore);
00648     cpl_vector_delete(sens_mini);
00649 
00650     return cpl_error_get_code();
00651 
00652 }
00653 
00654 
00655 /*----------------------------------------------------------------------------*/
00666 /*----------------------------------------------------------------------------*/
00667 static cpl_error_code visir_spc_phot_qc(cpl_propertylist * self,
00668                                         const char * star_name,
00669                                         double exptime,
00670                                         double sens_mean,
00671                                         double sens_stdev)
00672 {
00673 
00674     bug_if (cpl_propertylist_append_double(self, "ESO QC EXPTIME", exptime));
00675     bug_if (cpl_propertylist_append_double(self, "ESO QC SENS MEAN", sens_mean));
00676     bug_if (cpl_propertylist_append_double(self, "ESO QC SENS STDEV",
00677                                            sens_stdev));
00678     /* FIXME: Verify that this is the correct behaviour */
00679     bug_if (cpl_propertylist_append_string(self, "ESO QC STARNAME",
00680                                            star_name ? star_name : ""));
00681 
00682     end_skip;
00683 
00684     return cpl_error_get_code();
00685 
00686 }

Generated on Fri Jul 3 11:15:23 2009 for VISIR Pipeline Reference Manual by  doxygen 1.5.8