visir_spectro.c

00001 /* $Id: visir_spectro.c,v 1.205 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.205 $
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 <float.h>
00039 #include <assert.h>
00040 
00041 #include <cpl.h>
00042 
00043 #include "irplib_framelist.h"
00044 #include "irplib_plot.h"
00045 
00046 #include "visir_utils.h"
00047 #include "visir_pfits.h"
00048 #include "visir_inputs.h"
00049 #include "visir_spectro.h"
00050 
00051 /*----------------------------------------------------------------------------*/
00057 /*----------------------------------------------------------------------------*/
00058 
00059 /*-----------------------------------------------------------------------------
00060                             Private Function Prototypes
00061  -----------------------------------------------------------------------------*/
00062 static cpl_bivector * visir_spc_extract(cpl_image *, cpl_propertylist *,
00063                                         cpl_image **, int);
00064 
00065 static cpl_bivector * visir_bivector_load_fits(const char *, const char*,
00066                                                const char*);
00067 static cpl_error_code visir_bivector_interpolate(cpl_bivector *,
00068                                                  const cpl_bivector *);
00069 
00070 static cpl_error_code visir_spc_emission(cpl_bivector *, const cpl_vector *,
00071                                          const cpl_bivector *,
00072                                          const cpl_bivector *,
00073                                          const cpl_vector *, double);
00074 
00075 static cpl_polynomial * visir_spc_phys_disp(int, double, visir_spc_resol, int);
00076 static cpl_error_code visir_vector_convolve_symm(cpl_vector *,
00077                                                  const cpl_vector *);
00078 static cpl_image * visir_spc_flip(const cpl_image *, double, visir_spc_resol);
00079 static cpl_error_code visir_spc_xcorr(cpl_vector *, cpl_bivector *,
00080                                       cpl_vector *, const cpl_vector *,
00081                                       const cpl_bivector *,
00082                                       const cpl_bivector *,
00083                                       const cpl_vector *, const cpl_polynomial *,
00084                                       double, int, double, double *, int *);
00085 
00086 static cpl_vector * cpl_spc_convolve_init(int, double, double, int);
00087 
00088 static cpl_error_code visir_spectro_qclist_wcal(cpl_propertylist *,
00089                                                 int, double, double,
00090                                                 const cpl_polynomial *,
00091                                                 const cpl_polynomial *);
00092 
00093 static cpl_error_code visir_spectro_qclist_obs(cpl_propertylist *,
00094                                                double, double);
00095 
00096 static const double N_upper = 13.4e-6; /* Upper limit of N-band */
00097 static const double whechelle = 35.8/2; /* Half the echelle width */
00098 
00099 #ifndef VISIR_XC_LEN
00100 #define VISIR_XC_LEN 50
00101 #endif
00102 #ifndef VISIR_XC_FLEN
00103 #define VISIR_XC_FLEN 3
00104 #endif
00105 #ifndef VISIR_XC_SUBSEARCH
00106 #define VISIR_XC_SUBSEARCH 100
00107 #endif
00108 
00111 /*-----------------------------------------------------------------------------
00112                                 Function code
00113  -----------------------------------------------------------------------------*/
00114 
00115 /*----------------------------------------------------------------------------*/
00130 /*----------------------------------------------------------------------------*/
00131 visir_spc_resol visir_spc_get_res_wl(const irplib_framelist * rawframes,
00132                                      double * pwlen, double * pslitw,
00133                                      double * ptemp, double * pfwhm)
00134 {
00135     cpl_errorstate cleanstate = cpl_errorstate_get();
00136      /* Avoid (false) uninit warning */
00137     visir_spc_resol    resol = VISIR_SPC_R_ERR;
00138     char               ptmp[IRPLIB_FITS_STRLEN];
00139     double             wl, spx;
00140     double             sl = 0.0; /* Avoid (false) uninit warning */
00141     cpl_boolean        need_temp = ptemp != NULL;
00142     int                n;
00143     int                i;
00144 
00145     /* Check entries */
00146     cpl_ensure(rawframes != NULL, CPL_ERROR_NULL_INPUT, VISIR_SPC_R_ERR);
00147     cpl_ensure(pwlen     != NULL, CPL_ERROR_NULL_INPUT, VISIR_SPC_R_ERR);
00148     cpl_ensure(pslitw    != NULL, CPL_ERROR_NULL_INPUT, VISIR_SPC_R_ERR);
00149     cpl_ensure(pfwhm     != NULL, CPL_ERROR_NULL_INPUT, VISIR_SPC_R_ERR);
00150 
00151     n = irplib_framelist_get_size(rawframes);
00152 
00153     cpl_ensure(n > 0, CPL_ERROR_DATA_NOT_FOUND, VISIR_SPC_R_ERR);
00154 
00155      /* Allow 1 nm difference */
00156     skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_DOUBLE_WLEN,
00157                                       CPL_TYPE_DOUBLE, CPL_TRUE, 1e-3));
00158 
00159      /* Allow 1 micron difference */
00160     skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_DOUBLE_PIXSPACE,
00161                                       CPL_TYPE_DOUBLE, CPL_TRUE, 1e-6));
00162 
00163     /* The actual value depends on the age of the file :-( */
00164     skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_DOUBLE_SLITWIDTH,
00165                                       CPL_TYPE_DOUBLE, CPL_FALSE, 0.0));
00166 
00167     skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_STRING_RESOL,
00168                                       CPL_TYPE_STRING, CPL_TRUE, 0.0));
00169 
00170     skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_STRING_SLITNAME,
00171                                       CPL_TYPE_STRING, CPL_TRUE, 0.0));
00172 
00173     for (i=0; i < n; i++) {
00174         const cpl_propertylist * plist;
00175         const char * filename =
00176             cpl_frame_get_filename(irplib_framelist_get_const(rawframes, i));
00177         const char * pfits;
00178         double             wl_tmp, sl_tmp, spx_tmp;
00179 
00180 
00181         cpl_ensure(!cpl_error_get_code(), CPL_ERROR_DATA_NOT_FOUND,
00182                       VISIR_SPC_R_ERR);
00183 
00184         cpl_ensure(filename != NULL, CPL_ERROR_DATA_NOT_FOUND,
00185                       VISIR_SPC_R_ERR);
00186 
00187         plist = irplib_framelist_get_propertylist_const(rawframes, i);
00188 
00189         cpl_ensure(plist != NULL, CPL_ERROR_DATA_NOT_FOUND, VISIR_SPC_R_ERR);
00190 
00191         wl_tmp = visir_pfits_get_wlen(plist); 
00192         sl_tmp = visir_pfits_get_slitwidth(plist);
00193         spx_tmp = visir_pfits_get_pixspace(plist);
00194         pfits = visir_pfits_get_resol(plist);
00195         
00196         cpl_ensure(!cpl_error_get_code(), CPL_ERROR_DATA_NOT_FOUND,
00197                       VISIR_SPC_R_ERR);
00198 
00199         if (i == 0) {
00200             
00201             visir_optmod ins_settings;
00202 
00203             sl = sl_tmp;
00204             spx = spx_tmp;
00205             wl = wl_tmp;
00206 
00207             /* Divide the slit width with the
00208                Spectral PFOV = 0.127 Arcseconds/pixel */
00209             /* FIXME: The Spectral PFOV may change with a new detector */
00210             *pslitw = sl / 0.127; /* Convert Slit width from Arcseconds to pixel */
00211 
00212             *pwlen = wl * 1e-6; /* Convert from micron to m */
00213 
00214             strncpy(ptmp, pfits, IRPLIB_FITS_STRLEN);
00215             ptmp[IRPLIB_FITS_STRLEN] = '\0';
00216 
00217             cpl_msg_info(cpl_func, "RESOL [LR|MR|HRS|HRG] and WLEN [m] (%d frames)"
00218                          ": %s %g", n, ptmp, *pwlen);
00219 
00220             if (spx <= 0) {
00221                 cpl_msg_error(cpl_func,"Pixel Spacing (%g) in %s is non-positive",
00222                               spx, filename);
00223                 cpl_ensure(0, CPL_ERROR_ILLEGAL_INPUT, VISIR_SPC_R_ERR);
00224             }
00225 
00226             if (*pslitw <= 0) {
00227                 cpl_msg_error(cpl_func,"Slit Width (%g) in %s is non-positive",
00228                               sl, filename);
00229                 cpl_ensure(0, CPL_ERROR_ILLEGAL_INPUT, VISIR_SPC_R_ERR);
00230             }
00231 
00232             cpl_msg_info(cpl_func, "Slit Width [pixel] and Pixel Spacing [m]: "
00233                          "%g %g", *pslitw, spx);
00234 
00235             if (!strcmp("LR", ptmp)) {
00236                 resol = VISIR_SPC_R_LR;
00237             } else if (!strcmp("MR", ptmp)) {
00238                 resol = VISIR_SPC_R_MR;
00239             } else if (!strcmp("HRS", ptmp)) {
00240                 resol = VISIR_SPC_R_HR;
00241             } else if (!strcmp("HRG", ptmp)) {
00242                 resol = VISIR_SPC_R_GHR;
00243             } else {
00244                 cpl_msg_error(cpl_func,"Unsupported resolution (%s) in %s",
00245                               ptmp, filename);
00246                 cpl_ensure(0, CPL_ERROR_UNSUPPORTED_MODE, VISIR_SPC_R_ERR);
00247             }
00248             if (visir_spc_optmod_init(resol, *pwlen, &ins_settings)) {
00249                 cpl_msg_error(cpl_func, "Resolution %s does not support "
00250                               "Central Wavelength [m]: %g", ptmp, *pwlen);
00251                 cpl_ensure(0, CPL_ERROR_INCOMPATIBLE_INPUT, VISIR_SPC_R_ERR);
00252             }
00253 
00254             cpl_msg_info(cpl_func, "The %s-Spectral Resolution at %gm: %g",
00255                          ptmp, *pwlen,
00256                          visir_spc_optmod_resolution(&ins_settings));
00257             cpl_msg_info(cpl_func, "The %s-Linear Dispersion at %gm [pixel/m]: "
00258                          "%g", ptmp, *pwlen,
00259                          visir_spc_optmod_dispersion(&ins_settings));
00260 
00261             *pfwhm  = *pwlen * visir_spc_optmod_dispersion(&ins_settings)
00262                 / visir_spc_optmod_resolution(&ins_settings);
00263 
00264             cpl_msg_info(cpl_func, "The %s-FWHM at %gm [pixel]: %g",
00265                          ptmp, *pwlen, *pfwhm);
00266         } else {
00267             if (fabs(sl-sl_tmp) > 1e-3) { /* Allow 1 micron difference */
00268                 cpl_msg_error(cpl_func, "Inconsistent slit width (%g <=>"
00269                               " %g) in %s (%d of %d)",
00270                               sl, sl_tmp, filename, i+1, n);
00271                 cpl_ensure(0, CPL_ERROR_INCOMPATIBLE_INPUT, VISIR_SPC_R_ERR);
00272             }
00273         }
00274         if (need_temp) {
00275             /* Temperature [Celcius] not yet found */
00276             const double temp = visir_pfits_get_temp(plist);
00277             if (cpl_error_get_code()) {
00278                 visir_error_reset("Could not get FITS key");
00279             } else if ((-20 < temp) && (temp < 60)) {
00280                 /* Only accept a non-extreme temperature */
00281                 need_temp = CPL_FALSE;
00282                 *ptemp = temp;
00283             }
00284         }
00285 
00286     }
00287 
00288     if (need_temp) {
00289         cpl_msg_warning(cpl_func, "No FITS-files specify the M1 temperature, "
00290                      "using default");
00291         *ptemp = 10; /* Default is 10 Celcius */
00292     }
00293 
00294 
00295     if (ptemp != NULL) {
00296         *ptemp += 273.15; /* Convert to Kelvin */
00297         cpl_msg_info(cpl_func, "The M1 temperature [Kelvin]: %g", *ptemp);
00298     }
00299 
00300     end_skip;
00301 
00302     return resol;
00303 
00304 }
00305 
00306 /*----------------------------------------------------------------------------*/
00327 /*----------------------------------------------------------------------------*/
00328 cpl_error_code visir_vector_resample(cpl_vector * self, 
00329                                      const cpl_vector * xbounds,
00330                                      const cpl_bivector * source)
00331 {
00332 
00333     const cpl_vector * xsource  = cpl_bivector_get_x_const(source);
00334     const cpl_vector * ysource  = cpl_bivector_get_y_const(source);
00335 
00336     const double     * pxsource = cpl_vector_get_data_const(xsource);
00337     const double     * pysource = cpl_vector_get_data_const(ysource);
00338     const double     * pxbounds = cpl_vector_get_data_const(xbounds);
00339 
00340 
00341     cpl_vector   * ybounds  = cpl_vector_new(cpl_vector_get_size(xbounds));
00342     cpl_bivector * boundary = cpl_bivector_wrap_vectors((cpl_vector*)xbounds,
00343                                                         ybounds);
00344     double       * pybounds = cpl_vector_get_data(ybounds);
00345 
00346     double       * pself  = cpl_vector_get_data(self);
00347     const int      npix     = cpl_vector_get_size(self);
00348     int i;
00349     int itt;
00350 
00351 
00352     cpl_ensure_code(cpl_bivector_get_size(boundary) == npix + 1,
00353                         CPL_ERROR_ILLEGAL_INPUT);
00354 
00355     skip_if (0);
00356 
00357     itt = cpl_vector_find(xsource, pxbounds[0]);
00358 
00359     skip_if (0);
00360 
00361     skip_if (visir_bivector_interpolate(boundary, source));
00362 
00363     /* At this point itt most likely points to element just below
00364        pxbounds[0] */
00365     while (pxsource[itt] < pxbounds[0]) itt++;
00366 
00367     for (i=0; i < npix; i++) {
00368 
00369         /* The i'th value is the weighted average of the two interpolated
00370            values at the boundaries and the source values in between */
00371 
00372         double xlow  = pxbounds[i];
00373         double x     = pxsource[itt];
00374 
00375         if (x > pxbounds[i+1]) x = pxbounds[i+1];
00376         /* Contribution from interpolated value at lower boundary */
00377         pself[i] = pybounds[i] * (x - xlow);
00378 
00379         /* Contribution from table values in between boundaries */
00380         while (pxsource[itt] < pxbounds[i+1]) {
00381             const double xprev = x;
00382             x = pxsource[itt+1];
00383             if (x > pxbounds[i+1]) x = pxbounds[i+1];
00384             pself[i] += pysource[itt] * (x - xlow);
00385             xlow = xprev;
00386             itt++;
00387         }
00388 
00389         /* Contribution from interpolated value at upper boundary */
00390         pself[i] += pybounds[i+1] * (pxbounds[i+1] - xlow);
00391 
00392         /* Compute average by dividing integral by length of sampling interval
00393            (the factor 2 comes from the contributions) */
00394         pself[i] /= 2 * (pxbounds[i+1] - pxbounds[i]);
00395 
00396     }
00397 
00398 
00399     end_skip;
00400 
00401     cpl_vector_delete(ybounds);
00402     cpl_bivector_unwrap_vectors(boundary);
00403 
00404     return cpl_error_get_code();
00405 }
00406 
00407 
00408 
00409 /*----------------------------------------------------------------------------*/
00433 /*----------------------------------------------------------------------------*/
00434 cpl_error_code visir_spc_extract_wcal(const cpl_image * combined,
00435                                       const cpl_image * hcycle,
00436                                       double wlen, double slitw,
00437                                       double temp, double fwhm,
00438                                       visir_spc_resol resol,
00439                                       int ioffset,
00440                                       const char * spc_cal_lines,
00441                                       const char * spc_cal_qeff,
00442                                       cpl_table ** pspc_table,
00443                                       cpl_image ** pweight2d,
00444                                       cpl_propertylist * qclist,
00445                                       int doplot)
00446 {
00447 
00448     /* Both spectrum and error */
00449     cpl_bivector  * spc_n_err = NULL;
00450     cpl_image     * flipped   = NULL;
00451     const int       npix = cpl_image_get_size_y(combined);
00452 
00453 
00454     cpl_ensure_code(pweight2d != NULL, CPL_ERROR_NULL_INPUT);
00455 
00456     *pweight2d = NULL;
00457 
00458     cpl_ensure_code(npix > 0, CPL_ERROR_ILLEGAL_INPUT);
00459     cpl_ensure_code(npix == cpl_image_get_size_y(hcycle),
00460                         CPL_ERROR_ILLEGAL_INPUT);
00461 
00462 
00463     skip_if (0);
00464 
00465     skip_if (visir_spc_wavecal(hcycle, qclist, wlen, slitw, temp, fwhm, resol,
00466                                ioffset, spc_cal_lines, spc_cal_qeff,
00467                                pspc_table, doplot));
00468 
00469     /* Convert the combined image */
00470     flipped = visir_spc_flip(combined, wlen, resol);
00471     skip_if (0);
00472 
00473     /* Extract spectrum with error from the combined image */
00474     /* FIXME: Move inside */
00475     spc_n_err = visir_spc_extract(flipped, qclist, pweight2d,
00476                                   doplot);
00477     skip_if (0);
00478 
00479     cpl_image_delete(flipped);
00480     flipped = NULL;
00481 
00482     skip_if (*pspc_table == NULL);
00483 
00484     skip_if (cpl_table_new_column(*pspc_table, "SPC_EXTRACTED", CPL_TYPE_DOUBLE));
00485     skip_if (cpl_table_new_column(*pspc_table, "SPC_ERROR", CPL_TYPE_DOUBLE));
00486 
00487     skip_if (cpl_table_set_column_unit(*pspc_table, "SPC_EXTRACTED", "ADU/s"));
00488     skip_if (cpl_table_set_column_unit(*pspc_table, "SPC_ERROR", "ADU/s"));
00489 
00490     skip_if (cpl_table_copy_data_double(*pspc_table, "SPC_EXTRACTED", 
00491                                         cpl_bivector_get_x_data(spc_n_err)));
00492     skip_if (cpl_table_copy_data_double(*pspc_table, "SPC_ERROR", 
00493                                         cpl_bivector_get_y_data(spc_n_err)));
00494 
00495     if (doplot) {
00496         visir_table_plot("set grid;set xlabel 'Wavelength [m]';",
00497                          "t 'Extracted Spectrum' w linespoints",
00498                          "", *pspc_table, "WLEN", "SPC_EXTRACTED");
00499         visir_table_plot("set grid;set xlabel 'Wavelength [m]';",
00500                          "t 'Error on Extracted Spectrum' w linespoints",
00501                          "", *pspc_table, "WLEN", "SPC_ERROR");
00502     }
00503 
00504     end_skip;
00505 
00506     cpl_image_delete(flipped);
00507     cpl_bivector_delete(spc_n_err);
00508 
00509     return cpl_error_get_code();
00510 }
00511 
00512 
00513 /*----------------------------------------------------------------------------*/
00538 /*----------------------------------------------------------------------------*/
00539 cpl_error_code visir_spc_wavecal(const cpl_image * hcycle,
00540                                  cpl_propertylist * qclist,
00541                                  double wlen, double slitw,
00542                                  double temp, double fwhm,
00543                                  visir_spc_resol resol,
00544                                  int ioffset,
00545                                  const char * linefile,
00546                                  const char * qefffile,
00547                                  cpl_table ** pspc_table, int doplot)
00548 {
00549 
00550     /* Dispersion relation from physical model */
00551     cpl_polynomial * phdisp = NULL;
00552     /* Dispersion relation corrected by cross-correlation */
00553     cpl_polynomial * xcdisp = NULL;
00554 
00555     cpl_bivector * emission = NULL;
00556     cpl_vector   * boundary = NULL;
00557 
00558     cpl_bivector * temiss = NULL;
00559     cpl_bivector * tqeff  = NULL;
00560 
00561     cpl_image    * corrected = NULL;
00562 
00563     cpl_image    * xc_image  = NULL;
00564     cpl_vector   * xc_vector = NULL;
00565     cpl_bivector * xc_subres = NULL;
00566 
00567     cpl_vector   * vsymm   = NULL;
00568 
00569     cpl_vector   * vxc       = NULL;
00570 
00571     cpl_vector   * xc_subresx;
00572     cpl_vector   * xc_subresy;
00573 
00574     const int      npix = cpl_image_get_size_y(hcycle);
00575     int            delta, bestdelta, rawdelta;
00576     double         subdelta;
00577     double         xc0;
00578     double         qcxc, qcsubdelta;
00579     double         hc_min;
00580     int            convohlen;
00581     int            xc_flen;
00582     const int      i0 = 0;
00583     const int      i1 = 1;
00584     int            i;
00585     int            minpos;
00586     double       * pemiss;
00587     cpl_vector   * xemiss;
00588 
00589 
00590     assert( VISIR_XC_LEN >=0 && VISIR_XC_FLEN >=0);
00591     assert( VISIR_XC_SUBSEARCH == 1 ||
00592            (VISIR_XC_SUBSEARCH  > 1 && (VISIR_XC_SUBSEARCH&1)) == 0);
00593 
00594     cpl_ensure_code(!cpl_error_get_code(), cpl_error_get_code());
00595     cpl_ensure_code(pspc_table, CPL_ERROR_NULL_INPUT);
00596     cpl_ensure_code(npix > 0,   CPL_ERROR_ILLEGAL_INPUT);
00597 
00598 
00599     /* Make sure the corrected image is of type double */
00600     corrected = cpl_image_cast(hcycle, CPL_TYPE_DOUBLE);
00601     skip_if (0);
00602 
00603     hc_min = cpl_image_get_min(corrected);
00604     skip_if (0);
00605     cpl_msg_info(cpl_func,"Half-cycle image [%d X %d] has minimum intensity: %g",
00606                  cpl_image_get_size_x(hcycle), npix, hc_min);
00607     if (hc_min < 0) {
00608         cpl_msg_warning(cpl_func, "Thresholding negative intensities in half-"
00609                         "cycle image: %g", hc_min);
00610         skip_if (cpl_image_threshold(corrected, 0.0, DBL_MAX, 0.0, DBL_MAX));
00611     } else if (hc_min > 0) {
00612         skip_if (cpl_image_subtract_scalar(corrected, hc_min));
00613     }      
00614 
00615     /* Average the spatial dimension - into a cpl_vector */
00616     xc_image = cpl_image_collapse_create(corrected, 1);
00617     skip_if (0);
00618     skip_if (cpl_image_divide_scalar(xc_image, npix));
00619 
00620     cpl_image_delete(corrected);
00621     corrected = NULL;
00622 
00623     /* The dispersion relation goes from the top of the image to the bottom */
00624     if (resol == VISIR_SPC_R_HR || resol == VISIR_SPC_R_GHR) {
00625         /* Flip (if A-side), nothing else */
00626         corrected = visir_spc_flip(xc_image, wlen, resol);
00627         skip_if (0);
00628 
00629         cpl_image_delete(xc_image);
00630         xc_image = corrected;
00631         corrected = NULL;
00632     } else {
00633         skip_if (cpl_image_flip(xc_image, 0));
00634     }
00635 
00636     xc_vector = cpl_vector_wrap(npix, cpl_image_get_data(xc_image));
00637     skip_if (0);
00638 
00639     emission = cpl_bivector_new(npix + 2 * VISIR_XC_LEN);
00640     skip_if (0);
00641 
00642     boundary = cpl_vector_new(npix + 2 * VISIR_XC_LEN + 1);
00643     skip_if (0);
00644 
00645     phdisp = visir_spc_phys_disp(npix, wlen, resol, ioffset);
00646     skip_if (0);
00647 
00648     cpl_msg_info(cpl_func, "Dispersion polynomial of physical model:"
00649                  " %gm + ipix * %gm/pixel [ipix = 1, 2, ..., %d]",
00650                  cpl_polynomial_get_coeff(phdisp, &i0),
00651                  cpl_polynomial_get_coeff(phdisp, &i1), npix);
00652 
00653     temiss = visir_bivector_load_fits(linefile, "Wavelength", "Emission");
00654     if (cpl_error_get_code()) {
00655         cpl_msg_error(cpl_func, "Could not load file with Emission Lines");
00656         skip_if (1);
00657     }
00658 
00659     tqeff  = visir_bivector_load_fits(qefffile, "Wavelength", "Efficiency");
00660     if (cpl_error_get_code()) {
00661         cpl_msg_error(cpl_func, "Could not load file with Quantum-Efficiencies");
00662         skip_if (1);
00663     }
00664 
00665     *pspc_table = cpl_table_new(npix);
00666     skip_if (0);
00667 
00668     skip_if (cpl_table_new_column(*pspc_table, "WLEN", CPL_TYPE_DOUBLE));
00669     skip_if (cpl_table_new_column(*pspc_table, "SPC_MODEL_PH", CPL_TYPE_DOUBLE));
00670     skip_if (cpl_table_new_column(*pspc_table, "SPC_MODEL_XC", CPL_TYPE_DOUBLE));
00671     skip_if (cpl_table_new_column(*pspc_table, "SPC_SKY", CPL_TYPE_DOUBLE));
00672 
00673     skip_if (cpl_table_set_column_unit(*pspc_table, "WLEN", "m"));
00674     skip_if (cpl_table_set_column_unit(*pspc_table, "SPC_MODEL_PH",
00675                                        "J*radian/m^3/s"));
00676     skip_if (cpl_table_set_column_unit(*pspc_table, "SPC_MODEL_XC",
00677                                        "J*radian/m^3/s"));
00678     skip_if (cpl_table_set_column_unit(*pspc_table, "SPC_SKY", "ADU/s"));
00679 
00680 
00681     vsymm = cpl_spc_convolve_init(npix, slitw, fwhm, doplot);
00682 
00683     skip_if (vsymm == NULL);
00684 
00685     convohlen = cpl_vector_get_size(vsymm);
00686 
00687     skip_if (convohlen < 1);
00688 
00689     xc_flen = convohlen-1 < VISIR_XC_FLEN ? VISIR_XC_FLEN
00690         : (convohlen-1 > VISIR_XC_LEN ? VISIR_XC_LEN : convohlen-1);
00691 
00692 
00693     /* Determine the (possibly large) initial pixel shift */
00694   
00695     xc_subres = cpl_bivector_new(VISIR_XC_SUBSEARCH);
00696     skip_if (0);
00697 
00698     xc_subresy = cpl_bivector_get_y(xc_subres);
00699     skip_if (0);
00700     xc_subresx = cpl_bivector_get_x(xc_subres);
00701     skip_if (0);
00702 
00703     /* Copy the dispersion relation */
00704     xcdisp = cpl_polynomial_new(1);
00705     skip_if (cpl_polynomial_copy(xcdisp, phdisp));
00706 
00707 
00708     vxc = cpl_vector_new(2 * VISIR_XC_LEN + 1);
00709     skip_if (visir_spc_xcorr(vxc, emission, boundary, xc_vector, temiss, tqeff,
00710                              vsymm, xcdisp, -VISIR_XC_LEN, VISIR_XC_LEN,
00711                              temp, &qcxc, &rawdelta));
00712 
00713     if (doplot > 0) {
00714         cpl_vector   * xaxis = cpl_vector_new(2 * VISIR_XC_LEN + 1);
00715         cpl_bivector * bivxc = cpl_bivector_wrap_vectors(xaxis, vxc);
00716 
00717         for (i=0; i < 2 * VISIR_XC_LEN + 1; i++)
00718             if (cpl_vector_set(xaxis, i, i-VISIR_XC_LEN)) break;
00719 
00720         if (!cpl_error_get_code())
00721             visir_bivector_plot("set grid;set xlabel 'Offset [pixel]';",
00722                                 "t 'Cross-correlation (coarse)'", "", bivxc);
00723         cpl_bivector_unwrap_vectors(bivxc);
00724         cpl_vector_delete(xaxis);
00725     }
00726 
00727     skip_if (cpl_vector_set_size(vxc, 2 * VISIR_XC_FLEN + 1));
00728 
00729     skip_if (cpl_vector_set(xc_subresx, VISIR_XC_SUBSEARCH/2, rawdelta));
00730     skip_if (cpl_vector_set(xc_subresy, VISIR_XC_SUBSEARCH/2, qcxc));
00731 
00732     qcsubdelta = rawdelta;
00733     bestdelta = 0;
00734 
00735     cpl_msg_debug(cpl_func, "xc (%d): %g", rawdelta, qcxc);
00736 
00737     /*  Dump the unshifted model spectrum to the table
00738         - The unshifted signal starts at index VISIR_XC_LEN */
00739     pemiss = cpl_bivector_get_y_data(emission) + VISIR_XC_LEN;
00740     skip_if (cpl_table_copy_data_double(*pspc_table, "SPC_MODEL_PH", pemiss));
00741 
00742     /* Apply the initial pixel shift */
00743 #if defined CPL_VERSION_CODE && CPL_VERSION_CODE > CPL_VERSION(4, 5, 0)
00744     skip_if (cpl_polynomial_shift_1d(xcdisp, 0, rawdelta));
00745 #else
00746     skip_if (cpl_polynomial_shift_1d(xcdisp, rawdelta));
00747 #endif
00748 
00749     /* emission & boundary can be made shorter, but npix+VISIR_XC_FLEN
00750        elements must be free of edge-convolution effects */
00751     cpl_bivector_delete(emission);
00752     emission = NULL;
00753     cpl_vector_delete(boundary);
00754     boundary = NULL;
00755 
00756     emission = cpl_bivector_new(npix + 2 * xc_flen);
00757     skip_if (0);
00758 
00759     boundary = cpl_vector_new(npix + 2 * xc_flen + 1);
00760     skip_if (0);
00761 
00762     /* subdelta search starts with an offset of minus a half pixel
00763        and is in the range [-0.5; 0.5 [ */
00764     minpos = 0;
00765     subdelta = VISIR_XC_SUBSEARCH == 1 ? 0 : -0.5;
00766     for (i = 0; i < VISIR_XC_SUBSEARCH; i++,
00767              subdelta += 1/(double)VISIR_XC_SUBSEARCH) {
00768         double xc;
00769 
00770         if (2*i == VISIR_XC_SUBSEARCH) continue; /* subdelta == 0 */
00771 
00772         skip_if (visir_spc_xcorr(vxc, emission, boundary, xc_vector, temiss,
00773                                  tqeff, vsymm, xcdisp, -xc_flen + subdelta,
00774                                  VISIR_XC_FLEN, temp, &xc, &delta));
00775 
00776         skip_if (cpl_vector_set(xc_subresx, i, rawdelta+delta+subdelta));
00777         skip_if (cpl_vector_set(xc_subresy, i, xc));
00778         if (rawdelta+delta+subdelta < cpl_vector_get(xc_subresx, minpos))
00779             minpos = i;
00780 
00781         cpl_msg_debug(cpl_func, "xc (%g): %g %g", rawdelta+delta+subdelta, xc,
00782                       qcxc);
00783 
00784         skip_if (0);
00785 
00786         if (xc <= qcxc) continue; /* FIXME: Reverse expression ?! */
00787 
00788         qcxc = xc;
00789         bestdelta = delta;
00790         qcsubdelta = delta + subdelta + rawdelta;
00791 
00792     }
00793 
00794     if (minpos > 0) {
00795         /* Move the minimum offset to the beginning of the bivector */
00796         /* Currently only needed for plotting */
00797         const size_t size1 = sizeof(double) * minpos;
00798         const size_t size2 = sizeof(double) * (VISIR_XC_SUBSEARCH-minpos);
00799         double * swap = cpl_malloc(size1);
00800         double * pdata;
00801 
00802         pdata = cpl_vector_get_data(xc_subresx);
00803         memcpy(swap, pdata, size1);
00804         memmove(pdata, pdata + minpos, size2);
00805         memcpy(pdata+(VISIR_XC_SUBSEARCH-minpos), swap, size1);
00806 
00807         pdata = cpl_vector_get_data(xc_subresy);
00808         memcpy(swap, pdata, size1);
00809         memmove(pdata, pdata + minpos, size2);
00810         memcpy(pdata+(VISIR_XC_SUBSEARCH-minpos), swap, size1);
00811 
00812         cpl_free(swap);
00813     }
00814 
00815     cpl_vector_delete(boundary);
00816     boundary = NULL;
00817     cpl_bivector_delete(emission);
00818     emission = NULL;
00819 
00820     skip_if (0);
00821 
00822     if (fabs(qcsubdelta) >= VISIR_XC_LEN) {
00823         cpl_msg_warning(cpl_func, "Cross-correlation (%g pixel shift): %g",
00824                         qcsubdelta, qcxc);
00825     } else {
00826         cpl_msg_info(cpl_func,"XC pixel-shift: %d + %d + %g", rawdelta, bestdelta,
00827                      qcsubdelta - rawdelta - bestdelta);
00828         cpl_msg_info(cpl_func,"Cross-correlation (%g pixel shift): %g",
00829                      qcsubdelta, qcxc);
00830         assert( bestdelta <   VISIR_XC_LEN);
00831         assert( bestdelta >  -VISIR_XC_LEN);
00832     }
00833 
00834     if (qcxc <= 0) {
00835         /* Absolutely no cross-correlation */
00836         cpl_msg_error(cpl_func, "Atmospheric and Model Spectra have non-"
00837                       "positive cross-correlation (%g pixel shift): %g", 
00838                       qcsubdelta, qcxc);
00839         visir_error_set(CPL_ERROR_DATA_NOT_FOUND);
00840         skip_if(1);
00841     }
00842 
00843     /* Apply the sub-pixel precision shift - ignore the initial shift */
00844 #if defined CPL_VERSION_CODE && CPL_VERSION_CODE > CPL_VERSION(4, 5, 0)
00845     skip_if (cpl_polynomial_shift_1d(xcdisp, 0, qcsubdelta - rawdelta));
00846 #else
00847     skip_if (cpl_polynomial_shift_1d(xcdisp, qcsubdelta - rawdelta));
00848 #endif
00849 
00850     cpl_msg_info(cpl_func, "Dispersion polynomial from cross-correlation: "
00851                  "%gm + ipix * %gm/pixel [ipix = 1, 2, ..., %d]",
00852                  cpl_polynomial_get_coeff(xcdisp, &i0),
00853                  cpl_polynomial_get_coeff(xcdisp, &i1), npix);
00854 
00855     cpl_msg_info(cpl_func, "New Central Wavelength [m]: %g",
00856                  cpl_polynomial_eval_1d(xcdisp, 0.5*npix+0.5, NULL));
00857 
00858     /* Generate the new wavelengths based on the cross-correlation shift */
00859     emission = cpl_bivector_new(npix);
00860     xemiss = cpl_bivector_get_x(emission);
00861     skip_if (cpl_vector_fill_polynomial(xemiss, xcdisp, 1, 1));
00862 
00863     /* If the spectrum goes into N-band the sky spectrum may have variable
00864        atmospheric features, that are not in the model used for the model
00865        spectrum. This can cause the wavelength calibration to yield completely
00866        results */
00867     if (cpl_vector_get(xemiss,0) < N_upper &&
00868         N_upper < cpl_vector_get(xemiss,cpl_vector_get_size(xemiss)-1))
00869         cpl_msg_warning(cpl_func, "Spectrum goes above N-band (%gm). Wavelength "
00870                         "Calibration may be entirely inaccurate", N_upper);
00871 
00872     skip_if (cpl_table_copy_data_double(*pspc_table, "WLEN",
00873                                         cpl_bivector_get_x_data(emission)));
00874 
00875     /* - and the corresponding pixel boundaries */
00876     boundary = cpl_vector_new(npix + 1);
00877     skip_if (0);
00878     skip_if (cpl_vector_fill_polynomial(boundary, xcdisp, 0.5, 1));
00879 
00880     /* Get the emission at those wavelengths */
00881     skip_if (visir_spc_emission(emission, boundary, temiss, tqeff, vsymm,
00882                                 temp));
00883 
00884     skip_if (cpl_table_copy_data_double(*pspc_table, "SPC_MODEL_XC", 
00885                                         cpl_bivector_get_y_data(emission)));
00886 
00887     skip_if (cpl_table_copy_data_double(*pspc_table, "SPC_SKY", 
00888                                         cpl_vector_get_data(xc_vector)));
00889 
00890     /* The spectrum generated with xcdisp should have the maximum
00891        cross-correlation at zero offset */
00892     skip_if (cpl_vector_set_size(vxc, 1));
00893 
00894     delta = cpl_vector_correlate(vxc, cpl_bivector_get_y(emission),
00895                                  xc_vector);
00896     skip_if (delta < 0);
00897 
00898     xc0 = qcxc - cpl_vector_get(vxc, delta);
00899     cpl_vector_delete(vxc);
00900     vxc = NULL;
00901 
00902 #if 0
00903     /* FIXME: This check is broken with new concolution scheme */
00904     /* FIX ME: Why npix squared ? */
00905     /* The imperfect convolution at the spectral ends causes a warning here 
00906        when threshold is: 10 * npix * npix * DBL_EPSILON */
00907 
00908     if (delta || npix * fabs(xc0) > 25 * sigma)
00909         cpl_msg_warning(cpl_func, "Cross-correlation inconsistency(%d): %g",
00910                         delta, xc0);
00911 #endif
00912 
00913     if (doplot) {
00914         cpl_bivector * plot = cpl_bivector_wrap_vectors(xemiss,xc_vector);
00915 
00916         visir_bivector_plot("set grid;set xlabel 'Offset [pixel]';", "t 'Cross-"
00917                             "correlation (fine)' w linespoints", "", xc_subres);
00918 
00919         visir_bivector_plot("set grid;set xlabel 'Wavelength [m]';", "t 'Spec"
00920                             "trum from Half-cycle' w linespoints", "", plot);
00921         cpl_bivector_unwrap_vectors(plot);
00922 
00923         visir_bivector_plot("set grid;set xlabel 'Wavelength [m]';",
00924                              "t 'Shifted Model Spectrum' w linespoints",
00925                              "", emission);
00926 
00927         /* The unshifted model spectrum */
00928         visir_table_plot("set grid;set xlabel 'Wavelength [m]';",
00929                           "t 'Model Spectrum' w linespoints",
00930                           "", *pspc_table, "WLEN", "SPC_MODEL_PH");
00931 
00932     }
00933 
00934     /* Get the emissivity (range 0 to 1) for the calibrated wavelengths */
00935     skip_if (visir_vector_resample(cpl_bivector_get_y(emission),
00936                                       boundary, temiss));
00937 
00938     skip_if (cpl_table_new_column(*pspc_table, "SPC_EMISSIVITY",
00939                                   CPL_TYPE_DOUBLE));
00940 
00941     skip_if (cpl_table_copy_data_double(*pspc_table, "SPC_EMISSIVITY", 
00942                                         cpl_bivector_get_y_data(emission)));
00943 
00944     cpl_vector_delete(boundary);
00945     boundary = NULL;
00946 
00947     bug_if(visir_spectro_qclist_wcal(qclist, npix, qcxc, qcsubdelta,
00948                                      phdisp, xcdisp));
00949 
00950     if (doplot) {
00951 
00952         visir_bivector_plot("set grid;set xlabel 'Wavelength [m]';",
00953                              "t 'Atmospheric Emissivity' w linespoints",
00954                              "", emission);
00955 
00956         /* Create an model spectrum of twice the npix length */
00957         cpl_bivector_delete(emission);
00958         emission = cpl_bivector_new(2 * npix);
00959 
00960         boundary = cpl_vector_new(2 * npix + 1);
00961 
00962         cpl_vector_fill_polynomial(cpl_bivector_get_x(emission),
00963                                    phdisp, -0.5*npix, 1);
00964         cpl_vector_fill_polynomial(boundary, phdisp, -0.5*(npix+1), 1);
00965 
00966         /* Get the emission at those wavelengths */
00967         visir_spc_emission(emission, boundary, temiss, tqeff, vsymm, temp);
00968         cpl_vector_delete(boundary);
00969         boundary = NULL;
00970 
00971         visir_bivector_plot("set grid;set xlabel 'Wavelength [m]';",
00972                              "t 'Extended Model Spectrum' w linespoints",
00973                              "", emission);
00974 
00975     }
00976 
00977     end_skip;
00978 
00979     cpl_polynomial_delete(phdisp);
00980     cpl_polynomial_delete(xcdisp);
00981     cpl_image_delete(xc_image);
00982     cpl_vector_delete(vsymm);
00983     cpl_image_delete(corrected);
00984     cpl_bivector_delete(temiss);
00985     cpl_bivector_delete(tqeff);
00986     cpl_vector_delete(boundary);
00987     cpl_bivector_delete(emission);
00988     cpl_vector_unwrap(xc_vector);
00989     cpl_bivector_delete(xc_subres);
00990     cpl_vector_delete(vxc);
00991 
00992     return cpl_error_get_code();
00993 }
00994 
00995 
00996 /*----------------------------------------------------------------------------*/
01012 /*----------------------------------------------------------------------------*/
01013 cpl_error_code visir_spc_echelle_limit(int * pcol1, int * pcol2, double wlen,
01014                                        int ioffset, int icolmin, int icolmax)
01015 {
01016 
01017     visir_optmod ins_settings;
01018     double echpos;
01019     double wleni;   /* The central wavelength at order offset ioffset */
01020     int order;
01021     int error;
01022 
01023 
01024     cpl_ensure_code(wlen > 0,              CPL_ERROR_ILLEGAL_INPUT);
01025     cpl_ensure_code(pcol1,                 CPL_ERROR_NULL_INPUT);
01026     cpl_ensure_code(pcol2,                 CPL_ERROR_NULL_INPUT);
01027     cpl_ensure_code(icolmin > 0,           CPL_ERROR_ILLEGAL_INPUT);
01028     cpl_ensure_code(icolmax >= icolmin,    CPL_ERROR_ILLEGAL_INPUT);
01029     /* There are up to 5 spectra in the imaage */
01030     cpl_ensure_code(ioffset >= -4,         CPL_ERROR_ILLEGAL_INPUT);
01031     cpl_ensure_code(ioffset <=  4,         CPL_ERROR_ILLEGAL_INPUT);
01032 
01033     error = visir_spc_optmod_init(VISIR_SPC_R_GHR, wlen, &ins_settings);
01034     if (error) {
01035         cpl_msg_error(cpl_func, "HRG Optical model initialization (%p) failed: %d "
01036                       "(%g)", (void*)&ins_settings, error, wlen);
01037         cpl_ensure_code(0, CPL_ERROR_ILLEGAL_INPUT);
01038     }
01039     order = ioffset + visir_spc_optmod_get_echelle_order(&ins_settings);
01040 
01041     /* There are 18 echelle orders */
01042     cpl_ensure_code(order >   0,           CPL_ERROR_ILLEGAL_INPUT);
01043     cpl_ensure_code(order <= 18,           CPL_ERROR_ILLEGAL_INPUT);
01044 
01045     wleni = visir_spc_optmod_echelle(&ins_settings, wlen, ioffset  );
01046 
01047     echpos = visir_spc_optmod_cross_dispersion(&ins_settings, wleni);
01048     if (echpos <= whechelle || echpos >= icolmax-whechelle) {
01049         cpl_msg_error(cpl_func, "Echelle (%d) location out of range [%d;%d]: %g",
01050                       order, icolmin, icolmax, echpos);
01051         cpl_ensure_code(0, CPL_ERROR_DATA_NOT_FOUND);
01052     }
01053 
01054     *pcol1 = ceil(echpos - whechelle); /* Round up */
01055     *pcol2 = echpos + whechelle; /* Round down */
01056 
01057     if (*pcol1 < icolmin) *pcol1 = icolmin;
01058     if (*pcol2 > icolmax) *pcol2 = icolmax;
01059 
01060     cpl_msg_info(cpl_func, "Echelle order %d at col %g [%d; %d]", order, echpos,
01061                  *pcol1, *pcol2);
01062 
01063     return cpl_error_get_code();
01064 
01065 }
01066 
01067 /*----------------------------------------------------------------------------*/
01080 /*----------------------------------------------------------------------------*/
01081 cpl_image * visir_spc_column_extract(const cpl_image * self, int icol1,
01082                                      int icol2, int doplot)
01083 {
01084 
01085     cpl_image  * band    = NULL;
01086     cpl_image  * spatial = NULL;
01087     const int nrow = cpl_image_get_size_y(self);
01088     const int ncol = cpl_image_get_size_x(self);
01089 
01090     cpl_ensure(self != NULL,   CPL_ERROR_NULL_INPUT,    NULL);
01091     cpl_ensure(icol1 > 0,      CPL_ERROR_ILLEGAL_INPUT, NULL);
01092     cpl_ensure(icol2 >= icol1, CPL_ERROR_ILLEGAL_INPUT, NULL);
01093 
01094     cpl_ensure(ncol >= icol2,  CPL_ERROR_ILLEGAL_INPUT, NULL);
01095 
01096     band = cpl_image_extract(self, icol1, 1, icol2, nrow);
01097     skip_if (0);
01098 
01099     if (doplot > 0) {
01100         visir_image_plot("", "t 'The full-width image'", "", self);
01101 
01102         if (doplot > 1) {
01103             /* Average the spectral dimension */
01104             spatial = cpl_image_collapse_create(self, 0);
01105             skip_if (0);
01106             skip_if (cpl_image_divide_scalar(spatial, nrow));
01107 
01108             visir_image_row_plot("set grid;", "t 'Spectral direction "
01109                                  "collapsed' w linespoints", "",
01110                                  spatial, 1, 1, 1);
01111         }
01112     }
01113 
01114     end_skip;
01115 
01116     cpl_image_delete(spatial);
01117     if (cpl_error_get_code() && band != NULL) {
01118         cpl_image_delete(band);
01119         band = NULL;
01120     }
01121 
01122     return band;
01123 
01124 }
01125 
01126 
01127 /*----------------------------------------------------------------------------*/
01140 /*----------------------------------------------------------------------------*/
01141 cpl_error_code visir_spectro_qc(cpl_propertylist * qclist,
01142                                 cpl_propertylist * paflist,
01143                                 cpl_boolean        drop_wcs,
01144                                 const irplib_framelist * rawframes,
01145                                 const char * regcopy,
01146                                 const char * regcopypaf)
01147 {
01148 
01149     const cpl_propertylist * reflist
01150         = irplib_framelist_get_propertylist_const(rawframes, 0);
01151 
01152     bug_if (0);
01153 
01154     bug_if (visir_qc_append_capa(qclist, rawframes));
01155 
01156     if (regcopy != NULL)
01157         bug_if (cpl_propertylist_copy_property_regexp(qclist, reflist,
01158                                                       regcopy, 0));
01159 
01160     if (regcopypaf != NULL)
01161         bug_if (cpl_propertylist_copy_property_regexp(paflist, reflist,
01162                                                       regcopypaf, 0));
01163 
01164     bug_if (cpl_propertylist_append(paflist, qclist));
01165 
01166     if (drop_wcs) {
01167         cpl_propertylist * pcopy = cpl_propertylist_new();
01168         const cpl_error_code error
01169             = cpl_propertylist_copy_property_regexp(pcopy, reflist, "^("
01170                                                     IRPLIB_PFITS_WCS_REGEXP
01171                                                     ")$", 0);
01172         if (!error && cpl_propertylist_get_size(pcopy) > 0) {
01173             cpl_msg_warning(cpl_func, "Combined image will have no WCS "
01174                             "coordinates");
01175         }
01176         cpl_propertylist_delete(pcopy);
01177         bug_if(0);
01178     } else {
01179         bug_if(cpl_propertylist_copy_property_regexp(qclist, reflist, "^("
01180                                                      IRPLIB_PFITS_WCS_REGEXP
01181                                                      ")$", 0));
01182     }
01183 
01184     end_skip;
01185 
01186     return cpl_error_get_code();
01187 
01188 }
01189 
01190 
01194 /*----------------------------------------------------------------------------*/
01206 /*----------------------------------------------------------------------------*/
01207 static cpl_error_code visir_spectro_qclist_wcal(cpl_propertylist * self,
01208                                                 int npix, double xc,
01209                                                 double subdelta,
01210                                                 const cpl_polynomial * phdisp,
01211                                                 const cpl_polynomial * xcdisp)
01212 {
01213 
01214     const int    phdegree = cpl_polynomial_get_degree(phdisp);
01215     const int    xcdegree = cpl_polynomial_get_degree(xcdisp);
01216 
01217     const double phdisp1  = cpl_polynomial_get_coeff(phdisp, &phdegree);
01218     const double phdisp0  = cpl_polynomial_eval_1d(phdisp, 1.0, NULL);
01219 
01220     const double xcdisp1  = cpl_polynomial_get_coeff(xcdisp, &xcdegree);
01221     const double xcdisp0  = cpl_polynomial_eval_1d(xcdisp, 1.0, NULL);
01222 
01223     const double xcwlen   = cpl_polynomial_eval_1d(xcdisp, 0.5*(double)npix+0.5,
01224                                                    NULL);
01225 
01226 
01227     bug_if (0);
01228     bug_if (phdegree != 1);
01229     bug_if (xcdegree != 1);
01230 
01231     bug_if (cpl_propertylist_append_double(self, "ESO QC XC",       xc));
01232     bug_if (cpl_propertylist_append_double(self, "ESO QC XCSHIFT",  subdelta));
01233 
01234     bug_if (cpl_propertylist_append_int(self,    "ESO QC PHDEGREE", phdegree));
01235     bug_if (cpl_propertylist_append_double(self, "ESO QC PHDISPX1", phdisp1));
01236     bug_if (cpl_propertylist_append_double(self, "ESO QC PHDISPX0", phdisp0));
01237 
01238     bug_if (cpl_propertylist_append_double(self, "ESO QC XCWLEN",   xcwlen));
01239 
01240     bug_if (cpl_propertylist_append_int(self,    "ESO QC XCDEGREE", xcdegree));
01241     bug_if (cpl_propertylist_append_double(self, "ESO QC XCDISPX1", xcdisp1));
01242     bug_if (cpl_propertylist_append_double(self, "ESO QC XCDISPX0", xcdisp0));
01243 
01244     end_skip;
01245 
01246     return cpl_error_get_code();
01247 
01248 }
01249 
01250 
01251 
01252 /*----------------------------------------------------------------------------*/
01264 /*----------------------------------------------------------------------------*/
01265 static cpl_error_code visir_spectro_qclist_obs(cpl_propertylist * self,
01266                                                double xfwhm, double xcentro)
01267 {
01268 
01269 
01270     bug_if (0);
01271 
01272     bug_if (cpl_propertylist_append_double(self, "ESO QC XFWHM",    xfwhm));
01273     bug_if (cpl_propertylist_append_double(self, "ESO QC XCENTROI", xcentro));
01274 
01275     end_skip;
01276 
01277     return cpl_error_get_code();
01278 
01279 }
01280 
01281 
01282 /*----------------------------------------------------------------------------*/
01308 /*----------------------------------------------------------------------------*/
01309 static cpl_error_code visir_bivector_interpolate(cpl_bivector * out,
01310                                                  const cpl_bivector * in)
01311 {
01312     const cpl_error_code err = CPL_ERROR_ILLEGAL_INPUT;
01313 
01314     int m, n;
01315 
01316     const double * xref;
01317     const double * yref;
01318     double * xout;
01319     double * yout;
01320 
01321     /* Initialize to avoid unjustified compiler warning */
01322     double grad = 0.0;
01323     double y00  = 0.0;
01324     /* Start interpolation from below */
01325     int iabove = 0;
01326     int ibelow = 0;  /* Avoid (false) uninit warning */
01327     int i;
01328 
01329 
01330     cpl_ensure_code(out,   CPL_ERROR_NULL_INPUT);
01331     cpl_ensure_code(in,    CPL_ERROR_NULL_INPUT);
01332 
01333     m = cpl_bivector_get_size(in);
01334     n = cpl_bivector_get_size(out);
01335 
01336     cpl_ensure_code(m > 1, err);
01337     cpl_ensure_code(n > 0, err);
01338 
01339     xref = cpl_bivector_get_x_data_const(in);
01340     yref = cpl_bivector_get_y_data_const(in);
01341     xout = cpl_bivector_get_x_data(out);
01342     yout = cpl_bivector_get_y_data(out);
01343 
01344     assert( xref);
01345     assert( yref);
01346     assert( xout);
01347     assert( yout);
01348 
01349     /* Verify that extrapolation is not necessary */
01350     cpl_ensure_code(xref[0  ] <= xout[0  ], err);
01351     cpl_ensure_code(xout[0  ] <  xout[n-1], err);
01352     cpl_ensure_code(xout[n-1] <= xref[m-1], err);
01353 
01354     for (i = 0; i < n; i++) {
01355         /* When possible reuse reference function abscissa points */
01356         if (xout[i] > xref[iabove] || i == 0) {
01357             /* No, need new points */
01358             while (xout[i] > xref[++iabove]);
01359             ibelow = iabove - 1;
01360 
01361             /* Verify that reference abscissa points are valid */
01362             cpl_ensure_code(xref[iabove] > xref[ibelow], err);
01363 
01364             grad = (yref[iabove] - yref[ibelow])
01365                  / (xref[iabove] - xref[ibelow]);
01366 
01367             y00   = yref[ibelow] - grad * xref[ibelow];
01368         } else
01369             /* Interpolation point may not be smaller than
01370                the lower reference point */
01371             cpl_ensure_code(xout[i] >= xref[ibelow], err);
01372 
01373         yout[i] = y00 + grad * xout[i];
01374 
01375     }
01376 
01377     return CPL_ERROR_NONE;
01378 }
01379 
01380 /*----------------------------------------------------------------------------*/
01392 /*----------------------------------------------------------------------------*/
01393 static cpl_error_code visir_vector_convolve_symm(cpl_vector * self,
01394                                                  const cpl_vector * vsymm)
01395 {
01396 
01397     const int      npix = cpl_vector_get_size(self);
01398     const int      ihwidth = cpl_vector_get_size(vsymm) - 1;
01399     cpl_vector   * raw     = cpl_vector_duplicate(self);
01400     double       * pself= cpl_vector_get_data(self);
01401     double       * praw    = cpl_vector_get_data(raw);
01402     const double * psymm  = cpl_vector_get_data_const(vsymm);
01403 
01404     int i, j;
01405 
01406 
01407     skip_if (0);
01408 
01409     /* The convolution does not support this */
01410     skip_if (ihwidth >= npix);
01411 
01412     /* Convolve with the symmetric function */
01413     for (i = 0; i < ihwidth; i++) {
01414         pself[i] = praw[i] * psymm[0];
01415         for (j = 1; j <= ihwidth; j++) {
01416             const int k = i-j < 0 ? 0 : i-j;
01417             pself[i] += (praw[k]+praw[i+j]) * psymm[j];
01418         }
01419 
01420     }
01421 
01422     for (i = ihwidth; i < npix-ihwidth; i++) {
01423         pself[i] = praw[i] * psymm[0];
01424         for (j = 1; j <= ihwidth; j++)
01425             pself[i] += (praw[i-j]+praw[i+j]) * psymm[j];
01426 
01427     }
01428     for (i = npix-ihwidth; i < npix; i++) {
01429         pself[i] = praw[i] * psymm[0];
01430         for (j = 1; j <= ihwidth; j++) {
01431             const int k = i+j > npix-1 ? npix - 1 : i+j;
01432             pself[i] += (praw[k]+praw[i-j]) * psymm[j];
01433         }
01434 
01435     }
01436 
01437     end_skip;
01438 
01439     cpl_vector_delete(raw);
01440 
01441     return cpl_error_get_code();
01442 }
01443 
01444 /*----------------------------------------------------------------------------*/
01463 /*----------------------------------------------------------------------------*/
01464 static cpl_image * visir_spc_flip(const cpl_image * image, double wlen,
01465                               visir_spc_resol resol)
01466 {
01467     cpl_image  * flipped = cpl_image_cast(image, CPL_TYPE_DOUBLE);
01468     visir_optmod ins_settings;
01469 
01470 
01471     skip_if (0);
01472 
01473     if ((resol == VISIR_SPC_R_HR || resol == VISIR_SPC_R_GHR) &&
01474         visir_spc_optmod_init(resol, wlen, &ins_settings)) {
01475         visir_error_set(CPL_ERROR_ILLEGAL_INPUT);
01476         skip_if (1);
01477     }
01478 
01479     /* The dispersion relation goes from the top of the image to the bottom
01480        - except using the B-side (in high resolution) */
01481     if ((resol != VISIR_SPC_R_HR && resol != VISIR_SPC_R_GHR) ||
01482         visir_spc_optmod_side_is_A(&ins_settings) > 0) {
01483 
01484         cpl_msg_info(cpl_func, "Flipping image");
01485 
01486         skip_if (cpl_image_flip(flipped, 0));
01487     }
01488 
01489     end_skip;
01490 
01491     if (cpl_error_get_code() && flipped) {
01492         cpl_image_delete(flipped);
01493         flipped = NULL;
01494     }
01495 
01496     return flipped;
01497 
01498 }
01499 
01500 /*----------------------------------------------------------------------------*/
01515 /*----------------------------------------------------------------------------*/
01516 static cpl_polynomial * visir_spc_phys_disp(int npix, double wlen,
01517                                             visir_spc_resol resol, int ioffset)
01518 {
01519 
01520     cpl_polynomial * phdisp = NULL;
01521     visir_optmod     ins_settings;
01522 
01523     double dwl;
01524     double wlen0;
01525     double wlen1;
01526     double disp;
01527     const int i1 = 1;
01528     const int i0 = 0;
01529 
01530 
01531     cpl_ensure(resol,    CPL_ERROR_ILLEGAL_INPUT, NULL);
01532     cpl_ensure(wlen > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
01533     cpl_ensure(npix > 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
01534 
01535     /* Initialize instrument-specific settings
01536         - the resolution is not needed hereafter
01537        visir_spc_optmod_init() does itself not use the CPL-error system
01538           because it is also used in a non-CPL scope */
01539 
01540     cpl_ensure(!visir_spc_optmod_init(resol, wlen, &ins_settings),
01541                CPL_ERROR_ILLEGAL_INPUT, NULL);
01542 
01543     /* Get wavelength range (and corresponding central-wavelength)
01544        visir_spc_optmod_wlen() does not use the CPL-error system
01545          because it is also used in a non-CPL scope */
01546     dwl = visir_spc_optmod_wlen(&ins_settings, &wlen0, &wlen1);
01547 
01548     cpl_ensure(dwl >= 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
01549 
01550     /* Central-wavelength residual on Scan-Angle determination */
01551     dwl -= wlen;
01552     /* Warn if the residual exceeds twice the machine-precision */
01553     if (fabs(dwl) > 2*wlen*DBL_EPSILON) cpl_msg_warning(cpl_func, "Too large res"
01554         "idual in Scan-Angle determination [meps]: %g", dwl/DBL_EPSILON/wlen);
01555 
01556     if ((resol == VISIR_SPC_R_HR || resol == VISIR_SPC_R_GHR) &&
01557         !visir_spc_optmod_side_is_A(&ins_settings)) {
01558         const double swap = wlen1;
01559         wlen1 = wlen0;
01560         wlen0 = swap;
01561     }
01562     cpl_ensure(wlen1 > wlen0, CPL_ERROR_ILLEGAL_INPUT, NULL);
01563 
01564     /* Construct the 1st degree dispersion relation
01565        based on the physical model */
01566     phdisp = cpl_polynomial_new(1);
01567 
01568     /* The dispersion */
01569     disp = (wlen1-wlen0)/(npix-1);
01570 
01571     skip_if (0);
01572 
01573     skip_if (cpl_polynomial_set_coeff(phdisp, &i1, disp));
01574 
01575     skip_if (cpl_polynomial_set_coeff(phdisp, &i0, wlen0-disp));
01576 
01577     if ((resol == VISIR_SPC_R_HR || resol == VISIR_SPC_R_GHR) &&
01578         !visir_spc_optmod_side_is_A(&ins_settings)) {
01579         cpl_msg_info(cpl_func,"HR B-side WLMin, WLMax, Disp: %g %g %g", wlen0,
01580                      wlen1, cpl_polynomial_get_coeff(phdisp, &i1));
01581     } else {
01582         cpl_msg_info(cpl_func,"WLMin, WLMax, Disp: %g %g %g", wlen0, wlen1,
01583                      cpl_polynomial_get_coeff(phdisp, &i1));
01584     }
01585 
01586     if (resol == VISIR_SPC_R_GHR && ioffset != 0) {
01587         /* Another HRG Echelle order is requested
01588            - shift the 1st degree polynomial */
01589         const double dispi = visir_spc_optmod_echelle(&ins_settings,
01590                                 cpl_polynomial_get_coeff(phdisp, &i1), ioffset);
01591         const double wlen0i= visir_spc_optmod_echelle(&ins_settings,
01592                                 cpl_polynomial_get_coeff(phdisp, &i0), ioffset);
01593 
01594         skip_if (cpl_polynomial_set_coeff(phdisp, &i1, dispi));
01595 
01596         skip_if (cpl_polynomial_set_coeff(phdisp, &i0, wlen0i));
01597 
01598         cpl_msg_info(cpl_func, "WLc relative error(%d): %g", ioffset,
01599                      (wlen0i - cpl_polynomial_eval_1d(phdisp, 1, NULL))/wlen0i);
01600     }
01601 
01602 
01603     end_skip;
01604 
01605     if (cpl_error_get_code() && phdisp != NULL) {
01606         cpl_polynomial_delete(phdisp);
01607         phdisp = NULL;
01608     }
01609 
01610     return phdisp;
01611 
01612 }
01613 
01614 
01615 /*----------------------------------------------------------------------------*/
01639 /*----------------------------------------------------------------------------*/
01640 static cpl_error_code visir_spc_xcorr(cpl_vector * vxc,
01641                                       cpl_bivector * emission,
01642                                       cpl_vector   * boundary,
01643                                       const cpl_vector * xc_vector,
01644                                       const cpl_bivector * temiss,
01645                                       const cpl_bivector * tqeff,
01646                                       const cpl_vector   * vsymm,
01647                                       const cpl_polynomial * xcdisp,
01648                                       double firstpix,
01649                                       int half_search,
01650                                       double temp,
01651                                       double * pxc,
01652                                       int    * pdelta)
01653 {
01654 
01655     cpl_ensure_code(emission, CPL_ERROR_NULL_INPUT);
01656     cpl_ensure_code(boundary, CPL_ERROR_NULL_INPUT);
01657     cpl_ensure_code(xc_vector,CPL_ERROR_NULL_INPUT);
01658     cpl_ensure_code(temiss,   CPL_ERROR_NULL_INPUT);
01659     cpl_ensure_code(tqeff,    CPL_ERROR_NULL_INPUT);
01660     cpl_ensure_code(vsymm,  CPL_ERROR_NULL_INPUT);
01661     cpl_ensure_code(xcdisp,   CPL_ERROR_NULL_INPUT);
01662     cpl_ensure_code(pxc,      CPL_ERROR_NULL_INPUT);
01663     cpl_ensure_code(pdelta,   CPL_ERROR_NULL_INPUT);
01664 
01665 
01666     /* Compute the wavelengths of the spectrum
01667        according to the physical model */
01668     skip_if (cpl_vector_fill_polynomial(cpl_bivector_get_x(emission),
01669                                         xcdisp, firstpix+1, 1));
01670     skip_if (cpl_vector_fill_polynomial(boundary, xcdisp,
01671                                         firstpix+0.5, 1));
01672 
01673     /* Get the emission at those wavelengths */
01674     skip_if (visir_spc_emission(emission, boundary, temiss, tqeff, vsymm,
01675                                 temp));
01676 
01677     *pdelta = cpl_vector_correlate(vxc, cpl_bivector_get_y(emission),
01678                                    xc_vector);
01679     skip_if (*pdelta < 0);
01680 
01681     *pxc = cpl_vector_get(vxc, *pdelta);
01682 
01683     skip_if (0);
01684 
01685     *pdelta -= half_search;
01686 
01687     end_skip;
01688 
01689     return cpl_error_get_code();
01690 
01691 }
01692 
01693 
01694 /*----------------------------------------------------------------------------*/
01707 /*----------------------------------------------------------------------------*/
01708 
01709 static cpl_bivector * visir_bivector_load_fits(const char * file,
01710                                                const char * labelx,
01711                                                const char * labely)
01712 {
01713 
01714     cpl_bivector * result = NULL;
01715     cpl_table    * table  = NULL;
01716     double       * prow;
01717     int            nlines;
01718 
01719 
01720     skip_if (0);
01721 
01722     table = cpl_table_load(file, 1, 0);
01723     if (cpl_error_get_code()) {
01724         cpl_msg_error(cpl_func, "Could not load FITS table from file: %s",
01725                       file ? file : "<NULL>");
01726         skip_if (1);
01727     }
01728 
01729     nlines = cpl_table_get_nrow(table);
01730     skip_if (0);
01731 
01732     prow = cpl_table_get_data_double(table, labelx);
01733     skip_if (0);
01734 
01735     result = cpl_bivector_new(nlines);
01736     skip_if (0);
01737 
01738     skip_if (!memcpy(cpl_bivector_get_x_data(result), prow,
01739                      nlines * sizeof(double)));
01740 
01741     prow = cpl_table_get_data_double(table, labely);
01742     skip_if (0);
01743 
01744     skip_if (!memcpy(cpl_bivector_get_y_data(result), prow,
01745                      nlines * sizeof(double)));
01746 
01747     cpl_msg_info(cpl_func, "Read %d rows from %s [%g;%g]",
01748                  nlines, file,
01749                  cpl_vector_get(cpl_bivector_get_x(result), 0),
01750                  cpl_vector_get(cpl_bivector_get_x(result), nlines-1));
01751 
01752     end_skip;
01753 
01754     cpl_table_delete(table);
01755 
01756     if (result && cpl_error_get_code()) {
01757         cpl_bivector_delete(result);
01758         result = NULL;
01759     }
01760 
01761     return result;
01762 
01763 }
01764 
01765 
01766 /*----------------------------------------------------------------------------*/
01793 /*----------------------------------------------------------------------------*/
01794 static cpl_error_code visir_spc_emission(cpl_bivector       * emission,
01795                                          const cpl_vector   * boundary,
01796                                          const cpl_bivector * temiss,
01797                                          const cpl_bivector * tqeff,
01798                                          const cpl_vector   * vsymm,
01799                                          double temp)
01800 {
01801     cpl_bivector * tqeffi   = NULL;
01802     cpl_vector   * planck   = NULL;
01803     const int      npix = cpl_bivector_get_size(emission);
01804 
01805 
01806     cpl_ensure_code(emission, CPL_ERROR_NULL_INPUT);
01807     cpl_ensure_code(boundary, CPL_ERROR_NULL_INPUT);
01808     cpl_ensure_code(temiss,   CPL_ERROR_NULL_INPUT);
01809     cpl_ensure_code(tqeff,    CPL_ERROR_NULL_INPUT);
01810 
01811     /* npix is currently 256 */
01812     cpl_ensure_code(npix > 1, CPL_ERROR_ILLEGAL_INPUT);
01813 
01814     cpl_ensure_code(cpl_vector_get_size(boundary) == npix + 1,
01815                         CPL_ERROR_ILLEGAL_INPUT);
01816 
01817 
01818     skip_if (0);
01819 
01820     planck = cpl_vector_new(npix);
01821     skip_if (0);
01822 
01823     /* The atmospheric emission is assumed to be equivalent to that of
01824        a Black Body at 253 K */
01825     cpl_photom_fill_blackbody(planck, CPL_UNIT_ENERGYRADIANCE,
01826                               cpl_bivector_get_x(emission),
01827                               CPL_UNIT_LENGTH, 253);
01828 
01829     skip_if (visir_vector_resample(cpl_bivector_get_y(emission),
01830                                       boundary, temiss));
01831 
01832     /* Convolve to reflect the instrument resolution */
01833     skip_if (visir_vector_convolve_symm(cpl_bivector_get_y(emission),
01834                                             vsymm));
01835 
01836     skip_if (cpl_vector_multiply(cpl_bivector_get_y(emission), planck));
01837 
01838     /* The telescope emission is assumed to be equivalent to that of
01839        a Black Body */
01840     cpl_photom_fill_blackbody(planck, CPL_UNIT_ENERGYRADIANCE,
01841                               cpl_bivector_get_x(emission),
01842                               CPL_UNIT_LENGTH, temp);
01843 
01844     /* The telescope emissivity is assumed to be uniform at 0.12 */
01845     skip_if (cpl_vector_multiply_scalar(planck, 0.12));
01846 
01847     /* Add the telescope emission to the atmospheric */
01848     skip_if (cpl_vector_add(cpl_bivector_get_y(emission), planck));
01849 
01850     /* Multiply by the detector quantum efficiency */
01851     tqeffi = cpl_bivector_duplicate(emission);
01852     skip_if (visir_bivector_interpolate(tqeffi, tqeff));
01853 
01854     skip_if (cpl_vector_multiply(cpl_bivector_get_y(emission),
01855                                  cpl_bivector_get_y(tqeffi)));
01856 
01857     end_skip;
01858 
01859     cpl_bivector_delete(tqeffi);
01860     cpl_vector_delete(planck);
01861 
01862     return cpl_error_get_code();
01863 }
01864 
01865 
01866 /*----------------------------------------------------------------------------*/
01889 /*----------------------------------------------------------------------------*/
01890 static cpl_vector * cpl_spc_convolve_init(int maxlen, double slitw,
01891                                           double fwhm, int doplot)
01892 {
01893 
01894     const double sigma  = fwhm * CPL_MATH_SIG_FWHM;
01895     const int ihtophat  = (int)slitw/2;
01896     const int gausshlen = 1 + 5 * sigma + ihtophat < maxlen/2
01897         ? 1 + 5 * sigma + ihtophat : maxlen/2 - 1;
01898     /* convolen must be at least twice the gausshlen */
01899     const int convolen  = 1 + 10 * sigma + 8*ihtophat;
01900     cpl_vector * self = cpl_vector_new(gausshlen);
01901     cpl_vector * tophat = cpl_vector_new(convolen);
01902     int i;
01903 
01904     /* Easiest way to fill with a Gaussian is via a CPL image */
01905     cpl_image  * iself = cpl_image_wrap_double(gausshlen, 1,
01906                                                cpl_vector_get_data(self));
01907 
01908 
01909     skip_if (0);
01910 
01911     skip_if( slitw <= 0.0);
01912     skip_if( fwhm  <= 0.0);
01913     skip_if( convolen < 2 * gausshlen); /* This would indicate a bug */
01914 
01915     /* Place the top point of the Gaussian on left-most pixel */
01916     skip_if (cpl_image_fill_gaussian(iself, 1.0, 1.0, CPL_MATH_SQRT2PI,
01917                                      sigma, 1.0));
01918 
01919     if (doplot > 2) visir_vector_plot("set grid;", "t 'Right Half of Gaussian' "
01920                                       "w linespoints", "", self);
01921     
01922     /* The number of non-zero elements is 1+2*ihtophat */
01923     skip_if( cpl_vector_fill(tophat, 0.0));
01924 
01925     for (i = convolen/2-ihtophat; i < 1+convolen/2+ihtophat; i++)
01926         skip_if (cpl_vector_set(tophat, i, 1.0/(1.0+2.0*ihtophat)));
01927 
01928     /* Convolve the Top-hat with the Gaussian */
01929     skip_if (visir_vector_convolve_symm(tophat, self));
01930 
01931     if (doplot > 2) visir_vector_plot("set grid;","t 'Full Width Convolution' "
01932                                       "w linespoints", "", tophat);
01933     
01934     /* Overwrite the Gaussian with the Right Half of the convolution of the
01935        Top-hat + Gausssian */
01936 #if 1
01937     memcpy(cpl_vector_get_data(self),
01938            cpl_vector_get_data(tophat) + convolen/2,
01939            sizeof(double)*gausshlen);
01940 #else
01941     /* Equivalent, but slower */
01942     for (i = 0 ; i < gausshlen; i++)
01943         skip_if (cpl_vector_set(self, i, cpl_vector_get(tophat,
01944                                                           i + convolen/2)));
01945 #endif
01946 
01947     skip_if (0);
01948 
01949     cpl_msg_info(cpl_func, "Convolving Model Spectrum, Gauss-sigma=%g, "
01950                  "Tophat-width=%d, Truncation-Error=%g with width=%d", sigma,
01951                  1+2*ihtophat,
01952                  cpl_vector_get(self,gausshlen-1)/cpl_vector_get(self,0),
01953                  2*gausshlen-1);
01954 
01955     if (doplot > 1) visir_vector_plot("set grid;","t 'Right Half of Convolution"
01956                                       "' w linespoints", "", self);
01957 
01958     end_skip;
01959 
01960     cpl_vector_delete(tophat);
01961     cpl_image_unwrap(iself);
01962 
01963     if (cpl_error_get_code()) {
01964         cpl_vector_delete(self);
01965         self = NULL;
01966     }
01967 
01968     return self;
01969 
01970 }
01971 
01972 /*----------------------------------------------------------------------------*/
01989 /*----------------------------------------------------------------------------*/
01990 static cpl_bivector * visir_spc_extract(cpl_image * flipped,
01991                                         cpl_propertylist * qclist,
01992                                         cpl_image ** pweight2d,
01993                                         int doplot)
01994 {
01995     const int       ncol    = cpl_image_get_size_x(flipped);
01996     const int       npix    = cpl_image_get_size_y(flipped);
01997 
01998     cpl_bivector * result   = NULL;
01999     cpl_vector   * spectrum = NULL;
02000     cpl_vector   * error    = NULL;
02001     cpl_vector   * col      = NULL;
02002 
02003     cpl_image  * spatial  = NULL;
02004     cpl_image  * iweight  = NULL;
02005     cpl_vector * row      = NULL;
02006     cpl_image  * imrow    = NULL;
02007 
02008     double     * pweight  = NULL;
02009 
02010     cpl_apertures  * objects  = NULL;
02011     cpl_mask   * binary    = NULL;
02012     cpl_image  * locnoise  = NULL;
02013 
02014     double       xfwhm;   /* FWHM of brightest object */
02015     double       xcentro; /* X-Centroid of brightest object */
02016 
02017     int i, j;
02018     int is_rejected;
02019 
02020     const double sigma = 3; /* Assume signal at this level */
02021     double stdev2d, min, max, yfwhm;
02022     double weight_2norm;
02023     /* Position of the widest signal region */
02024     int ifwhm, jfwhm;
02025     int mspix;
02026     /* Low and High pixel of the widest signal-less region */
02027     int ilnoise, ihnoise;
02028     const int is_echelle = ncol <= 2 * (whechelle + 1);
02029 
02030 
02031     cpl_ensure(pweight2d != NULL, CPL_ERROR_NULL_INPUT, NULL);
02032 
02033     *pweight2d = NULL;
02034 
02035     skip_if (0);
02036 
02037     /* Compute spatial weights:
02038        mean-subtract each row and average + normalize */
02039 
02040     if (!is_echelle) {
02041         /* All but HR Grism has a negative signal equal to the positive
02042            i.e. the mean is zero */
02043         cpl_msg_info(cpl_func, "Combined image has mean: %g",
02044                      cpl_image_get_mean(flipped));
02045 
02046         col = cpl_vector_new(npix);
02047         skip_if (0);
02048 
02049         /* Subtract the mean from each row/wavelength */
02050         pweight = cpl_image_get_data(flipped);
02051         for (j=0; j < npix; j++, pweight += ncol) {
02052             double mean;
02053 
02054             imrow = cpl_image_wrap_double(1, ncol, pweight);
02055             skip_if (0);
02056 
02057             mean = cpl_image_get_mean(imrow);
02058             skip_if (0);
02059 
02060             skip_if (cpl_vector_set(col, j, mean));
02061 
02062             skip_if (cpl_image_subtract_scalar(imrow, mean));
02063 
02064             cpl_image_unwrap(imrow);
02065             imrow = NULL;
02066 
02067         }
02068 
02069         if (doplot > 1) visir_vector_plot("set grid;","t 'Estimated Background'"
02070                                           " w linespoints", "", col);
02071         cpl_vector_delete(col);
02072         col = NULL;
02073     }
02074 
02075     /* The st.dev. of the noise */
02076     stdev2d = visir_img_phot_sigma_clip(flipped)/sqrt(npix);
02077     skip_if (0);
02078 
02079     cpl_msg_info(cpl_func, "St.Dev. on noise in 2D-combined image: %g",
02080                  stdev2d);
02081 
02082     /* Average the spectral dimension */
02083     spatial = cpl_image_collapse_create(flipped, 0);
02084     skip_if (0);
02085     skip_if (cpl_image_divide_scalar(spatial, npix));
02086 
02087     skip_if (cpl_image_get_maxpos(spatial, &ifwhm, &jfwhm));
02088 
02089     if (doplot > 1) {
02090         visir_image_col_plot("","t 'Most intense column' w linespoints",
02091                              "", flipped, ifwhm, ifwhm, 1);
02092         visir_image_row_plot("set grid;", "t 'Combined image with "
02093                              "spectral direction collapsed' w linespoints",
02094                              "", spatial, 1, 1, 1);
02095     }
02096 
02097     max = cpl_image_get(spatial, ifwhm, 1, &is_rejected);
02098     if (is_rejected != 0 || max <= 0) {
02099         /* assert(is_rejected == 0); */
02100         cpl_msg_error(cpl_func, "Cannot compute FWHM on a collapsed spectrum "
02101                       "with a non-positive maximum: %g", max);
02102         visir_error_set(CPL_ERROR_DATA_NOT_FOUND);
02103         skip_if (1);
02104     }
02105 
02106     skip_if (cpl_image_get_fwhm(spatial, ifwhm, 1, &xfwhm,  &yfwhm));
02107 
02108     /* Remove noise from spatial */
02109     iweight = cpl_image_duplicate(spatial);
02110     skip_if (cpl_image_threshold(spatial, sigma*stdev2d,
02111                                  DBL_MAX, 0.0, DBL_MAX)
02112              || cpl_image_threshold(iweight, -FLT_MAX, -sigma*stdev2d,
02113                                     -FLT_MAX, 0.0));
02114 
02115     /* spatial has no negative values - find centroid */
02116     for (ilnoise = ifwhm; ilnoise > 0 &&
02117              cpl_image_get(spatial, ilnoise, 1, &is_rejected) > 0 &&
02118              !is_rejected; ilnoise--);
02119     skip_if (0);
02120     for (ihnoise = ifwhm; ihnoise <= ncol &&
02121              cpl_image_get(spatial, ihnoise, 1, &is_rejected) > 0 &&
02122              !is_rejected; ihnoise++);
02123     skip_if (0);
02124     /* There may be no negative weights at all */
02125     if (!ilnoise) ilnoise = 1;
02126     if (ihnoise > ncol) ihnoise = ncol;
02127 
02128     xcentro = cpl_image_get_centroid_x_window(spatial, ilnoise, 1, ihnoise, 1);
02129 
02130     cpl_msg_info(cpl_func, "Spatial FWHM(%d:%d:%d:%g): %g", ilnoise, ifwhm,
02131                  ihnoise, xcentro, xfwhm);
02132 
02133     skip_if (0);
02134     skip_if (cpl_image_add(spatial, iweight));
02135 
02136     skip_if (cpl_image_copy(iweight, spatial, 1, 1));
02137 
02138     /* Create weights that have an absolute sum of 1 - as an image */
02139     skip_if (cpl_image_normalise(iweight, CPL_NORM_ABSFLUX));
02140 
02141     if (doplot > 1) visir_image_row_plot("set grid;", "t 'Cleaned, normalized "
02142                                          "combined image with spectral direction"
02143                                          " averaged' w linespoints", "",
02144                                          iweight, 1, 1, 1);
02145 
02146     if (cpl_image_get_min(spatial) == 0 &&
02147         cpl_image_get_max(spatial) == 0) {
02148         visir_error_set(CPL_ERROR_DATA_NOT_FOUND);
02149         cpl_msg_error(cpl_func, "Spatial weights too noisy");
02150         skip_if (1);
02151     }
02152 
02153     weight_2norm = sqrt(cpl_image_get_sqflux(iweight));
02154 
02155     cpl_msg_info(cpl_func, "2-norm of weights: %g", weight_2norm);
02156 
02157 
02158 
02159     /* Determine st.dev. on noise at signal-less pixels */
02160     if (is_echelle) {
02161         int ileft = 5;
02162         int iright = ncol - 5;
02163         cpl_binary * pbin;
02164 
02165 
02166         if (ileft  > xcentro - xfwhm * 2)
02167             ileft  = xcentro - xfwhm * 2;
02168         if (iright < xcentro + xfwhm * 2)
02169             iright = xcentro + xfwhm * 2;
02170 
02171         cpl_msg_info(cpl_func, "HRG pixels of noise: [1 %d] [%d %d]", ileft,
02172                      iright, ncol);
02173 
02174         binary = cpl_mask_new(ncol, 1);
02175         skip_if (0);
02176 
02177         pbin = cpl_mask_get_data(binary);
02178         skip_if (0);
02179 
02180         for (i = 0; i < ncol; i++) pbin[i] = CPL_BINARY_0;
02181         for (i = 0; i < ileft; i++) pbin[i] = CPL_BINARY_1;
02182         for (i = iright; i < ncol; i++) pbin[i] = CPL_BINARY_1;
02183 
02184     } else {
02185         binary = cpl_mask_threshold_image_create(spatial, -sigma * stdev2d,
02186                                                  sigma * stdev2d);
02187     }
02188     skip_if (0);
02189 
02190     mspix = cpl_mask_count(binary);
02191     cpl_msg_info(cpl_func, "Pixels of noise(%g*%g): %d", stdev2d, sigma,
02192                  mspix);
02193     skip_if (0);
02194 
02195     if (mspix < 2) {
02196         /* No noise pixels found */
02197         cpl_msg_error(cpl_func, "Cannot estimate spectrum noise with just %d "
02198                       "pixels of noise", mspix);
02199         visir_error_set(CPL_ERROR_DATA_NOT_FOUND);
02200         skip_if (1);
02201     }
02202 
02203     locnoise = cpl_image_new_from_mask(binary);
02204     cpl_mask_delete(binary);
02205     binary = NULL;
02206 
02207     skip_if (0);
02208 
02209     error = cpl_vector_new(npix);
02210     skip_if (0);
02211 
02212 
02213     /* Compute for each wavelength the noise */
02214     for (j=0; j < npix; j++) {
02215 
02216         double npp, stdev1d;
02217 
02218 
02219         imrow = cpl_image_extract(flipped, 1, j+1, ncol, j+1);
02220 
02221         skip_if (0);
02222 
02223         objects = cpl_apertures_new_from_image(imrow, locnoise);
02224         cpl_image_delete(imrow);
02225         imrow = NULL;
02226              
02227         skip_if (0);
02228 
02229         stdev1d = cpl_apertures_get_stdev(objects, 1);
02230         cpl_apertures_delete(objects);
02231         objects = NULL;
02232 
02233         /* The noise per pixel is defined as the Standard Deviation
02234            on the noise (computed from the part of the signal that
02235            has no object signal) multiplied by the 2-norm of the
02236            noise-thresholded spatial weights */
02237 
02238         npp = weight_2norm * stdev1d;
02239 
02240         skip_if (cpl_vector_set(error, j, npp));
02241     }
02242 
02243     /* Spectrum noise computation done */
02244 
02245 
02246     /* Iterate through the spatial dimension - sum up the weighted column */
02247     for (i=1; i <= ncol; i++) {
02248         const double weight = cpl_image_get(iweight, i, 1, &is_rejected);
02249 
02250         skip_if (0);
02251         if (is_rejected) {
02252             /* This would require a whole column to be rejected */
02253             visir_error_set(CPL_ERROR_DATA_NOT_FOUND);
02254             skip_if (1);
02255         }
02256             
02257         /* The sigma-clipping may cause many columns to be zero */
02258         if (weight == 0) continue;
02259 
02260         col = cpl_vector_new_from_image_column(flipped, i); /* or medcorr */
02261         skip_if (0);
02262 
02263         skip_if (cpl_vector_multiply_scalar(col, weight));
02264 
02265         if (spectrum == NULL) {
02266             spectrum = col;
02267         } else {
02268             skip_if (cpl_vector_add(spectrum, col));
02269             cpl_vector_delete(col);
02270         }
02271         col = NULL;
02272     }
02273 
02274     /* assert( spectrum ); */
02275 
02276     min = cpl_vector_get_min(spectrum);
02277     if (min <0) cpl_msg_warning(cpl_func, "Extracted spectrum has negative "
02278                                 "intensity: %g", min);
02279 
02280     /* Create 2D-weight map by replicating the 1D-weights over the
02281        wavelengths */
02282 
02283     *pweight2d = cpl_image_new(ncol, npix, CPL_TYPE_DOUBLE);
02284 
02285     for (j=1; j <= npix; j++)
02286         skip_if (cpl_image_copy(*pweight2d, iweight, 1, j));
02287 
02288     if (doplot > 0) visir_image_plot("", "t 'The weight map'", "", *pweight2d);
02289 
02290     bug_if(visir_spectro_qclist_obs(qclist, xfwhm, xcentro));
02291 
02292     end_skip;
02293 
02294     cpl_image_delete(locnoise);
02295     cpl_mask_delete(binary);
02296     cpl_image_delete(spatial);
02297     cpl_apertures_delete(objects);
02298     cpl_vector_delete(col);
02299     cpl_vector_delete(row);
02300     cpl_image_delete(imrow);
02301     cpl_image_delete(iweight);
02302 
02303     if (cpl_error_get_code()) {
02304         cpl_vector_delete(spectrum);
02305         cpl_vector_delete(error);
02306     } else {
02307 
02308         result = cpl_bivector_wrap_vectors(spectrum, error);
02309 
02310         if (doplot > 2) visir_bivector_plot("", "t 'error versus spectrum'",
02311                                             "", result);
02312     }
02313 
02314     return result;
02315 }
02316 

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