irplib_wlxcorr.c

00001 /* $Id: irplib_wlxcorr.c,v 1.52 2009/10/21 14:50:20 llundin Exp $
00002  *
00003  * This file is part of the IRPLIB package
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., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00019  */
00020 
00021 /*
00022  * $Author: llundin $
00023  * $Date: 2009/10/21 14:50:20 $
00024  * $Revision: 1.52 $
00025  * $Name: HEAD $
00026  */
00027 
00028 #ifdef HAVE_CONFIG_H
00029 #include <config.h>
00030 #endif
00031 
00032 /*-----------------------------------------------------------------------------
00033                                 Includes
00034  -----------------------------------------------------------------------------*/
00035 
00036 #include <math.h>
00037 #include <string.h>
00038 
00039 #include <cpl.h>
00040 
00041 #include "irplib_wavecal_impl.h"
00042 
00043 #include "irplib_wlxcorr.h"
00044 
00045 /*----------------------------------------------------------------------------*/
00055 /*----------------------------------------------------------------------------*/
00056 
00057 /*-----------------------------------------------------------------------------
00058                            Defines
00059  -----------------------------------------------------------------------------*/
00060 
00061 #ifndef inline
00062 #define inline /* inline */
00063 #endif
00064 
00065 #define IRPLIB_MAX(A,B) ((A) > (B) ? (A) : (B))
00066 #define IRPLIB_MIN(A,B) ((A) < (B) ? (A) : (B))
00067 
00068 #define IRPLIB_PTR_SWAP(a,b)                                               \
00069     do { void * irplib_ptr_swap =(a);(a)=(b);(b)=irplib_ptr_swap; } while (0)
00070 
00071 /*-----------------------------------------------------------------------------
00072                            Private functions
00073  -----------------------------------------------------------------------------*/
00074 
00075 static void irplib_wlxcorr_estimate(cpl_vector *, cpl_vector *,
00076                                     const cpl_vector *,
00077                                     const cpl_bivector *,
00078                                     const cpl_vector *,
00079                                     const cpl_polynomial *,
00080                                     double, double);
00081 
00082 static int irplib_wlxcorr_signal_resample(cpl_vector *, const cpl_vector *, 
00083         const cpl_bivector *) ;
00084 static cpl_error_code cpl_vector_fill_lss_profile_symmetric(cpl_vector *,
00085                                                             double, double);
00086 static cpl_error_code irplib_wlcalib_fill_spectrum(cpl_vector *,
00087                                                 const cpl_bivector *,
00088                                                 const cpl_vector *,
00089                                                 const cpl_polynomial *, int);
00090 
00091 static cpl_boolean irplib_wlcalib_is_lines(const cpl_vector *,
00092                                         const cpl_polynomial *,
00093                                         int, double);
00094 
00098 /*----------------------------------------------------------------------------*/
00134 /*----------------------------------------------------------------------------*/
00135 cpl_polynomial * irplib_wlxcorr_best_poly(const cpl_vector     * spectrum,
00136                                           const cpl_bivector   * lines_catalog,
00137                                           int                    degree,
00138                                           const cpl_polynomial * guess_poly,
00139                                           const cpl_vector     * wl_error,
00140                                           int                    nsamples,
00141                                           double                 slitw,
00142                                           double                 fwhm,
00143                                           double               * xc,
00144                                           cpl_table           ** wlres,
00145                                           cpl_vector          ** xcorrs)
00146 {
00147     const int         spec_sz = cpl_vector_get_size(spectrum);
00148     const int         nfree   = cpl_vector_get_size(wl_error);
00149     int               ntests  = 1;
00150     cpl_vector      * model;
00151     cpl_vector      * vxc;
00152     cpl_vector      * init_pts_wl;
00153     cpl_matrix      * init_pts_x;
00154     cpl_vector      * pts_wl;
00155     cpl_vector      * vxcorrs;
00156     cpl_vector      * conv_kernel = NULL;    
00157     cpl_polynomial  * poly_sol;
00158     cpl_polynomial  * poly_candi;
00159     const double    * pwl_error = cpl_vector_get_data_const(wl_error); 
00160     const double    * dxc;
00161     const cpl_boolean symsamp = CPL_TRUE; /* init_pts_x is symmetric */
00162     const cpl_boolean is_lines
00163         = irplib_wlcalib_is_lines(cpl_bivector_get_x_const(lines_catalog),
00164                                guess_poly, spec_sz, 1.0);
00165     int               i;
00166 
00167     /* FIXME: Need mode parameter for catalogue type (lines <=> profile) */
00168 
00169     /* In case of failure */
00170     if (wlres  != NULL) *wlres  = NULL;
00171     if (xcorrs != NULL) *xcorrs = NULL;
00172 
00173     /* Useful for knowing if resampling is used */
00174     cpl_msg_debug(cpl_func, "Checking %d^%d dispersion polynomials (slitw=%g, "
00175                   "fwhm=%g) against %d-point observed spectrum with%s "
00176                   "catalog resampling", nsamples, nfree, slitw, fwhm, spec_sz,
00177                   is_lines ? "out" : "");
00178 
00179     cpl_ensure(xc            != NULL, CPL_ERROR_NULL_INPUT,    NULL);
00180     *xc = -1.0;
00181     cpl_ensure(spectrum      != NULL,  CPL_ERROR_NULL_INPUT,    NULL);
00182     cpl_ensure(lines_catalog != NULL, CPL_ERROR_NULL_INPUT,    NULL);
00183     cpl_ensure(guess_poly    != NULL, CPL_ERROR_NULL_INPUT,    NULL);
00184     cpl_ensure(wl_error      != NULL, CPL_ERROR_NULL_INPUT,    NULL);
00185     cpl_ensure(nfree         >= 2,    CPL_ERROR_ILLEGAL_INPUT, NULL);
00186     cpl_ensure(nsamples      >  0,    CPL_ERROR_ILLEGAL_INPUT, NULL);
00187     /* FIXME: degree is redundant */
00188     cpl_ensure(1 + degree   == nfree, CPL_ERROR_ILLEGAL_INPUT, NULL);
00189 
00190     cpl_ensure(cpl_polynomial_get_dimension(guess_poly) == 1,
00191                CPL_ERROR_ILLEGAL_INPUT, NULL);
00192 
00193     if (nsamples > 1) {
00194         /* Search place must consist of more than one point */
00195         /* FIXME: The bounds should probably not be negative */
00196         for (i = 0; i < nfree; i++) {
00197             if (pwl_error[i] != 0.0) break;
00198         }
00199         cpl_ensure(i < nfree, CPL_ERROR_ILLEGAL_INPUT, NULL);
00200     }
00201  
00202     if (!is_lines) {
00203         /* Create the convolution kernel */
00204         conv_kernel = irplib_wlxcorr_convolve_create_kernel(slitw, fwhm);
00205         cpl_ensure(conv_kernel   != NULL, CPL_ERROR_ILLEGAL_INPUT, NULL);
00206     }
00207 
00208     /* Create initial test points */
00209     init_pts_x  = cpl_matrix_new(1, nfree);
00210     init_pts_wl = cpl_vector_new(nfree);
00211     pts_wl      = cpl_vector_new(nfree);
00212     for (i = 0; i < nfree; i++) {
00213         const double xpos  = spec_sz * i / (double)degree;
00214         const double wlpos = cpl_polynomial_eval_1d(guess_poly, xpos, NULL)
00215             - 0.5 * pwl_error[i];
00216 
00217         cpl_matrix_set(init_pts_x, 0, i, xpos);
00218         cpl_vector_set(init_pts_wl,   i, wlpos);
00219 
00220         ntests *= nsamples; /* Count number of tests */
00221 
00222     }
00223 
00224     vxcorrs = xcorrs != NULL ? cpl_vector_new(ntests) : NULL;
00225 
00226     poly_sol   = cpl_polynomial_new(1);
00227     poly_candi = cpl_polynomial_new(1);
00228     model = cpl_vector_new(spec_sz);
00229     vxc = cpl_vector_new(1);
00230     dxc = cpl_vector_get_data_const(vxc);
00231    
00232     /* Create the polynomial candidates and estimate them */
00233     for (i=0; i < ntests; i++) {
00234         int    idiv = i;
00235         int    deg;
00236 
00237         /* Update wavelength at one anchor point - and reset wavelengths
00238            to their default for any anchor point(s) at higher wavelengths */
00239         for (deg = degree; deg >= 0; deg--, idiv /= nsamples) {
00240             const int imod = idiv % nsamples;
00241             const double wlpos = cpl_vector_get(init_pts_wl, deg)
00242                                + imod * pwl_error[deg] / nsamples;
00243 
00244             /* FIXME: If wlpos causes pts_wl to be non-increasing, the
00245                solution will be non-physical with no need for evaluation.
00246                (*xc could be set to -1 in this case). */
00247             cpl_vector_set(pts_wl, deg, wlpos);
00248 
00249             if (imod > 0) break;
00250         }
00251 
00252         /* Generate */
00253         cpl_polynomial_fit(poly_candi, init_pts_x, &symsamp, pts_wl,
00254                            NULL, CPL_FALSE, NULL, &degree);
00255         /* *** Estimate *** */
00256         irplib_wlxcorr_estimate(vxc, model, spectrum, lines_catalog,
00257                                 conv_kernel, poly_candi, slitw, fwhm);
00258         if (vxcorrs != NULL) cpl_vector_set(vxcorrs, i, *dxc);
00259         if (*dxc > *xc) {
00260             /* Found a better solution */
00261             *xc = *dxc;
00262             IRPLIB_PTR_SWAP(poly_sol, poly_candi);
00263         }
00264     }
00265 
00266     cpl_vector_delete(model);
00267     cpl_vector_delete(vxc);
00268     cpl_vector_delete(conv_kernel);
00269     cpl_vector_delete(pts_wl);
00270     cpl_matrix_delete(init_pts_x);
00271     cpl_vector_delete(init_pts_wl);
00272     cpl_polynomial_delete(poly_candi);
00273 
00274 #ifdef CPL_WLCALIB_FAIL_ON_CONSTANT
00275     /* FIXME: */
00276     if (cpl_polynomial_get_degree(poly_sol) == 0) {
00277         cpl_polynomial_delete(poly_sol);
00278         cpl_vector_delete(vxcorrs);
00279         *xc = 0.0;
00280         cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_OUTPUT,
00281                                     __FILE__, __LINE__, "Found a constant "
00282                                     "dispersion");
00283             cpl_errorstate_dump(prestate, CPL_FALSE, NULL);
00284         return NULL;
00285     }
00286 #endif
00287     
00288     if (wlres != NULL) {
00289         /* FIXME: A failure in the table creation is not considered a failure
00290            of the whole function call (although all outputs may be useless) */
00291 
00292         cpl_errorstate prestate = cpl_errorstate_get();
00293         /* Create the spc_table  */
00294         *wlres = irplib_wlxcorr_gen_spc_table(spectrum, lines_catalog, slitw,
00295                                               fwhm, guess_poly, poly_sol);
00296         if (*wlres == NULL) {
00297             cpl_polynomial_delete(poly_sol);
00298             cpl_vector_delete(vxcorrs);
00299             *xc = -1.0;
00300             cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_OUTPUT,
00301                                         __FILE__, __LINE__, "Cannot generate "
00302                                         "infos table");
00303             /* cpl_errorstate_dump(prestate, CPL_FALSE, NULL); */
00304             cpl_errorstate_set(prestate);
00305             return NULL;
00306         }
00307     } 
00308     
00309     if (xcorrs != NULL) {
00310         *xcorrs = vxcorrs;
00311     } else {
00312         /* assert(vxcorrs == NULL); */
00313     }
00314 
00315     return poly_sol;
00316 }
00317 
00318 /*----------------------------------------------------------------------------*/
00336 /*----------------------------------------------------------------------------*/
00337 cpl_table * irplib_wlxcorr_gen_spc_table(
00338         const cpl_vector        *   spectrum,
00339         const cpl_bivector      *   lines_catalog,
00340         double                      slitw,
00341         double                      fwhm,
00342         const cpl_polynomial    *   guess_poly,
00343         const cpl_polynomial    *   corr_poly)
00344 {
00345 
00346     cpl_vector      *   conv_kernel = NULL;
00347     cpl_bivector    *   gen_init ;
00348     cpl_bivector    *   gen_corr ;
00349     cpl_table       *   spc_table ;
00350     const double    *   pgen ;
00351     const double        xtrunc = 0.5 * slitw + 5.0 * fwhm * CPL_MATH_SIG_FWHM;
00352     const int           spec_sz = cpl_vector_get_size(spectrum);
00353     const cpl_boolean   guess_resamp
00354         = !irplib_wlcalib_is_lines(cpl_bivector_get_x_const(lines_catalog),
00355                                 guess_poly, spec_sz, 1.0);
00356     const cpl_boolean   corr_resamp
00357         = !irplib_wlcalib_is_lines(cpl_bivector_get_x_const(lines_catalog),
00358                                 corr_poly, spec_sz, 1.0);
00359     cpl_error_code      error;
00360 
00361     cpl_msg_debug(cpl_func, "Tabel for guess dispersion polynomial (slitw=%g, "
00362                   "fwhm=%g) with %d-point observed spectrum with%s catalog re"
00363                   "sampling", slitw, fwhm, spec_sz, guess_resamp ? "out" : "");
00364     cpl_msg_debug(cpl_func, "Tabel for corr. dispersion polynomial (slitw=%g, "
00365                   "fwhm=%g) with %d-point observed spectrum with%s catalog re"
00366                   "sampling", slitw, fwhm, spec_sz, corr_resamp ? "out" : "");
00367 
00368     /* Test inputs */
00369     cpl_ensure(spectrum, CPL_ERROR_NULL_INPUT, NULL) ;
00370     cpl_ensure(lines_catalog, CPL_ERROR_NULL_INPUT, NULL) ;
00371     cpl_ensure(guess_poly, CPL_ERROR_NULL_INPUT, NULL) ;
00372     cpl_ensure(corr_poly, CPL_ERROR_NULL_INPUT, NULL) ;
00373 
00374     /* Create the convolution kernel */
00375     if (guess_resamp || corr_resamp) {
00376         conv_kernel = irplib_wlxcorr_convolve_create_kernel(slitw, fwhm);
00377 
00378         if (conv_kernel == NULL) {
00379             cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
00380                                         __FILE__, __LINE__, "Cannot create "
00381                                         "convolution kernel") ;
00382             return NULL ;
00383         }
00384     }
00385 
00386     /* Get the emission at initial wavelengths */
00387     gen_init = cpl_bivector_new(spec_sz);
00388     if (guess_resamp) {
00389         error = irplib_wlcalib_fill_spectrum(cpl_bivector_get_y(gen_init),
00390                                           lines_catalog, conv_kernel,
00391                                           guess_poly, 0);
00392     } else {
00393         error = irplib_vector_fill_line_spectrum_model
00394             (cpl_bivector_get_y(gen_init), NULL, NULL,
00395              guess_poly, lines_catalog,
00396              slitw, fwhm, xtrunc, 0, CPL_FALSE, CPL_FALSE, NULL);
00397     }
00398 
00399     if (error || cpl_vector_fill_polynomial(cpl_bivector_get_x(gen_init),
00400                                             guess_poly, 1, 1)) {
00401         cpl_vector_delete(conv_kernel);
00402         cpl_bivector_delete(gen_init);
00403         cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
00404                                     __FILE__, __LINE__, "Cannot get the "
00405                                     "emission spectrum");
00406         return NULL;
00407     }
00408  
00409     /* Get the emission at corrected wavelengths */
00410     gen_corr = cpl_bivector_new(spec_sz);
00411     if (corr_resamp) {
00412         error = irplib_wlcalib_fill_spectrum(cpl_bivector_get_y(gen_corr),
00413                                           lines_catalog, conv_kernel,
00414                                           corr_poly, 0);
00415     } else {
00416         error = irplib_vector_fill_line_spectrum_model
00417             (cpl_bivector_get_y(gen_corr), NULL, NULL,
00418              corr_poly, lines_catalog,
00419              slitw, fwhm, xtrunc, 0, CPL_FALSE, CPL_FALSE, NULL);
00420     }
00421 
00422     if (error || cpl_vector_fill_polynomial(cpl_bivector_get_x(gen_corr),
00423                                             corr_poly, 1, 1)) {
00424         cpl_vector_delete(conv_kernel);
00425         cpl_bivector_delete(gen_init);
00426         cpl_bivector_delete(gen_corr) ;
00427         cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
00428                                     __FILE__, __LINE__, "Cannot get the "
00429                                     "emission spectrum");
00430         return NULL;
00431     }
00432     cpl_vector_delete(conv_kernel) ;
00433 
00434     /* Create the ouput table */
00435     spc_table = cpl_table_new(spec_sz);
00436     cpl_table_new_column(spc_table, IRPLIB_WLXCORR_COL_WAVELENGTH, 
00437             CPL_TYPE_DOUBLE);
00438     cpl_table_new_column(spc_table, IRPLIB_WLXCORR_COL_CAT_INIT, 
00439             CPL_TYPE_DOUBLE);
00440     cpl_table_new_column(spc_table, IRPLIB_WLXCORR_COL_CAT_FINAL, 
00441             CPL_TYPE_DOUBLE);
00442     cpl_table_new_column(spc_table, IRPLIB_WLXCORR_COL_OBS, CPL_TYPE_DOUBLE);
00443     
00444     /* Update table */
00445     pgen = cpl_bivector_get_x_data_const(gen_corr) ;
00446     cpl_table_copy_data_double(spc_table, IRPLIB_WLXCORR_COL_WAVELENGTH, pgen) ;
00447     pgen = cpl_bivector_get_y_data_const(gen_corr) ;
00448     cpl_table_copy_data_double(spc_table, IRPLIB_WLXCORR_COL_CAT_FINAL, pgen) ;
00449     pgen = cpl_vector_get_data_const(spectrum) ;
00450     cpl_table_copy_data_double(spc_table, IRPLIB_WLXCORR_COL_OBS, pgen) ;
00451     pgen = cpl_bivector_get_y_data_const(gen_init) ;
00452     cpl_table_copy_data_double(spc_table, IRPLIB_WLXCORR_COL_CAT_INIT, pgen);
00453     cpl_bivector_delete(gen_init);
00454     cpl_bivector_delete(gen_corr);
00455 
00456     return spc_table ;
00457 }
00458 
00459 /*----------------------------------------------------------------------------*/
00471 /*----------------------------------------------------------------------------*/
00472 cpl_bivector * irplib_wlxcorr_cat_extract(
00473         const cpl_bivector  *   lines_catalog,
00474         double                  wave_min,
00475         double                  wave_max)
00476 {
00477     const int           nlines = cpl_bivector_get_size(lines_catalog);
00478     int                 wave_min_id, wave_max_id ;
00479     cpl_vector       *  sub_cat_wl ;
00480     cpl_vector       *  sub_cat_int ;
00481     const cpl_vector *  xlines  = cpl_bivector_get_x_const(lines_catalog);
00482     const double     *  dxlines = cpl_vector_get_data_const(xlines);
00483 
00484     cpl_ensure(lines_catalog != NULL, CPL_ERROR_NULL_INPUT,    NULL);
00485 
00486     /* Find the 1st line */
00487     wave_min_id = cpl_vector_find(xlines, wave_min);
00488     /* The first line must be greater than (at least?) wave_min */
00489     if (dxlines[wave_min_id] <= wave_min) wave_min_id++;
00490 
00491     /* Find the last line */
00492     wave_max_id = cpl_vector_find(xlines, wave_max);
00493     /* The last line must be less than wave_max */
00494     if (dxlines[wave_max_id] >= wave_min) wave_max_id--;
00495 
00496     /* Checking the wavelength range at this point via the indices also
00497        verifies that they were not found using non-increasing wavelengths */
00498     cpl_ensure(wave_min_id <= wave_max_id, CPL_ERROR_ILLEGAL_INPUT, NULL);
00499 
00500     if (wave_min_id < 0 || wave_max_id == nlines) {
00501         cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
00502                                     __FILE__, __LINE__, "The %d-line catalogue "
00503                                     "has no lines in the range %g -> %g",
00504                                     nlines, wave_min, wave_max);
00505         return NULL ;
00506     }
00507 
00508     sub_cat_wl = cpl_vector_extract(xlines, wave_min_id, wave_max_id, 1);
00509     sub_cat_int = cpl_vector_extract(cpl_bivector_get_y_const(lines_catalog), 
00510                                      wave_min_id, wave_max_id, 1);
00511  
00512     return cpl_bivector_wrap_vectors(sub_cat_wl, sub_cat_int);
00513 }
00514 
00515 /*----------------------------------------------------------------------------*/
00532 /*----------------------------------------------------------------------------*/
00533 cpl_vector * irplib_wlxcorr_convolve_create_kernel(double  slitw,
00534                                                    double  fwhm)
00535 {
00536     const double  sigma  = fwhm * CPL_MATH_SIG_FWHM;
00537     const int     size   = 1 + (int)(5.0 * sigma + 0.5*slitw);
00538     cpl_vector  * kernel = cpl_vector_new(size);
00539 
00540 
00541     if (cpl_vector_fill_lss_profile_symmetric(kernel, slitw, fwhm)) {
00542         cpl_vector_delete(kernel);
00543         cpl_ensure(0, cpl_error_get_code(), NULL);
00544     }
00545 
00546     return kernel;
00547 }
00548 
00549 /*----------------------------------------------------------------------------*/
00562 /*----------------------------------------------------------------------------*/
00563 int irplib_wlxcorr_convolve(
00564         cpl_vector          *   smoothed,
00565         const cpl_vector    *   conv_kernel)
00566 {
00567     int             nsamples ;
00568     int             ihwidth ;
00569     cpl_vector  *   raw ;
00570     double      *   psmoothe ;
00571     double      *   praw ;
00572     const double*   psymm ;
00573     int             i, j ;
00574 
00575     /* Test entries */
00576     cpl_ensure(smoothed, CPL_ERROR_NULL_INPUT, -1) ;
00577     cpl_ensure(conv_kernel, CPL_ERROR_NULL_INPUT, -1) ;
00578 
00579     /* Initialise */
00580     nsamples = cpl_vector_get_size(smoothed) ;
00581     ihwidth = cpl_vector_get_size(conv_kernel) - 1 ;
00582     cpl_ensure(ihwidth<nsamples, CPL_ERROR_ILLEGAL_INPUT, -1) ;
00583     psymm = cpl_vector_get_data_const(conv_kernel) ;
00584     psmoothe = cpl_vector_get_data(smoothed) ;
00585     
00586     /* Create raw vector */
00587     raw = cpl_vector_duplicate(smoothed) ;
00588     praw = cpl_vector_get_data(raw) ;
00589 
00590     /* Convolve with the symmetric function */
00591     for (i=0 ; i<ihwidth ; i++) {
00592         psmoothe[i] = praw[i] * psymm[0];
00593         for (j=1 ; j <= ihwidth ; j++) {
00594             const int k = i-j < 0 ? 0 : i-j;
00595             psmoothe[i] += (praw[k]+praw[i+j]) * psymm[j];
00596         }
00597     }
00598 
00599     for (i=ihwidth ; i<nsamples-ihwidth ; i++) {
00600         psmoothe[i] = praw[i] * psymm[0];
00601         for (j=1 ; j<=ihwidth ; j++)
00602             psmoothe[i] += (praw[i-j]+praw[i+j]) * psymm[j];
00603     }
00604     for (i=nsamples-ihwidth ; i<nsamples ; i++) {
00605         psmoothe[i] = praw[i] * psymm[0];
00606         for (j=1 ; j<=ihwidth ; j++) {
00607             const int k = i+j > nsamples-1 ? nsamples - 1 : i+j;
00608             psmoothe[i] += (praw[k]+praw[i-j]) * psymm[j];
00609         }
00610     }
00611     cpl_vector_delete(raw) ;
00612     return 0 ;
00613 }
00614 
00615 /*----------------------------------------------------------------------------*/
00625 /*----------------------------------------------------------------------------*/
00626 int irplib_wlxcorr_plot_solution(
00627         const cpl_polynomial    *   init,
00628         const cpl_polynomial    *   comp,
00629         const cpl_polynomial    *   sol,
00630         int                         pix_start,
00631         int                         pix_stop) 
00632 {
00633     int                 nsamples, nplots ;
00634     cpl_vector      **  vectors ;
00635     cpl_bivector    *   bivector ;
00636     double              diff ;
00637     int                 i ;
00638     
00639     /* Test entries */
00640     if (init == NULL || comp == NULL) return -1 ;
00641     
00642     /* Initialise */
00643     nsamples = pix_stop - pix_start + 1 ;
00644     if (sol != NULL)    nplots = 3 ;
00645     else                nplots = 2 ;
00646     
00647     /* Create vectors */
00648     vectors = cpl_malloc((nplots+1)*sizeof(cpl_vector*)) ;
00649     for (i=0 ; i<nplots+1 ; i++) vectors[i] = cpl_vector_new(nsamples) ;
00650 
00651     /* First plot with the lambda/pixel relation */
00652     /* Fill vectors */
00653     for (i=0 ; i<nsamples ; i++) {
00654         cpl_vector_set(vectors[0], i, pix_start+i) ;
00655         cpl_vector_set(vectors[1], i, 
00656                 cpl_polynomial_eval_1d(init, (double)(pix_start+i), NULL)) ;
00657         cpl_vector_set(vectors[2], i, 
00658                 cpl_polynomial_eval_1d(comp, (double)(pix_start+i), NULL)) ;
00659         if (sol != NULL) 
00660             cpl_vector_set(vectors[3], i, 
00661                     cpl_polynomial_eval_1d(sol, (double)(pix_start+i), NULL)) ;
00662     }
00663 
00664     /* Plot */
00665     cpl_plot_vectors("set grid;set xlabel 'Position (pixels)';", 
00666         "t '1-Initial / 2-Computed / 3-Solution' w lines", 
00667         "", (const cpl_vector **)vectors, nplots+1);
00668 
00669     /* Free vectors */
00670     for (i=0 ; i<nplots+1 ; i++) cpl_vector_delete(vectors[i]) ;
00671     cpl_free(vectors) ;
00672 
00673     /* Allocate vectors */
00674     nplots -- ;
00675     vectors = cpl_malloc((nplots+1)*sizeof(cpl_vector*)) ;
00676     for (i=0 ; i<nplots+1 ; i++) vectors[i] = cpl_vector_new(nsamples) ;
00677     
00678     /* Second plot with the delta-lambda/pixel relation */
00679     /* Fill vectors */
00680     for (i=0 ; i<nsamples ; i++) {
00681         cpl_vector_set(vectors[0], i, pix_start+i) ;
00682         diff = cpl_polynomial_eval_1d(comp, (double)(pix_start+i), NULL) -
00683             cpl_polynomial_eval_1d(init, (double)(pix_start+i), NULL) ;
00684         cpl_vector_set(vectors[1], i, diff) ;
00685         if (sol != NULL) {
00686             diff = cpl_polynomial_eval_1d(sol, (double)(pix_start+i), NULL) -
00687                 cpl_polynomial_eval_1d(init, (double)(pix_start+i), NULL) ;
00688             cpl_vector_set(vectors[2], i, diff) ;
00689         }
00690     }
00691 
00692     /* Plot */
00693     if (sol == NULL) {
00694         bivector = cpl_bivector_wrap_vectors(vectors[0], vectors[1]) ;
00695         cpl_plot_bivector(
00696 "set grid;set xlabel 'Position (pixels)';set ylabel 'Wavelength difference';", 
00697             "t 'Computed-Initial wavelenth' w lines", "", bivector);
00698         cpl_bivector_unwrap_vectors(bivector) ;
00699     } else {
00700         cpl_plot_vectors("set grid;set xlabel 'Position (pixels)';", 
00701             "t '1-Computed - Initial / 2--Solution - Initial' w lines", 
00702             "", (const cpl_vector **)vectors, nplots+1);
00703     }
00704     
00705     /* Free vectors */
00706     for (i=0 ; i<nplots+1 ; i++) cpl_vector_delete(vectors[i]) ;
00707     cpl_free(vectors) ;
00708 
00709     /* Return */
00710     return 0 ;
00711 }
00712 
00713 /*----------------------------------------------------------------------------*/
00720 /*----------------------------------------------------------------------------*/
00721 int irplib_wlxcorr_plot_spc_table(
00722         const cpl_table     *   spc_table, 
00723         const char          *   title) 
00724 {
00725     char                title_loc[1024] ;
00726     cpl_vector      **  vectors ;
00727     cpl_vector      **  sub_vectors ;
00728     cpl_vector      *   tmp_vec ;
00729     int                 nsamples ;
00730     double              hsize_nm, max, mean1, mean3 ;
00731     int                 start_ind, stop_ind, nblines, hsize_pix ;
00732     int                 i, j ;
00733 
00734     /* Test entries */
00735     if (spc_table == NULL) return -1 ;
00736     
00737     /* Initialise */
00738     nsamples = cpl_table_get_nrow(spc_table) ;
00739     hsize_nm = 0.2 ;
00740     hsize_pix = 10 ;
00741     nblines = 0 ;
00742     sprintf(title_loc, 
00743         "t '%s - 1-Initial catalog/2-Corrected catalog/3-Observed' w lines",
00744         title) ;
00745     title_loc[1023] = (char)0 ;
00746     
00747     vectors = cpl_malloc(4*sizeof(cpl_vector*)) ;
00748     vectors[0] = cpl_vector_wrap(nsamples, 
00749             cpl_table_get_data_double((cpl_table*)spc_table,
00750                 IRPLIB_WLXCORR_COL_WAVELENGTH));
00751     vectors[1] = cpl_vector_wrap(nsamples, 
00752             cpl_table_get_data_double((cpl_table*)spc_table, 
00753                 IRPLIB_WLXCORR_COL_CAT_INIT));
00754     vectors[2] = cpl_vector_wrap(nsamples, 
00755             cpl_table_get_data_double((cpl_table*)spc_table, 
00756                 IRPLIB_WLXCORR_COL_CAT_FINAL));
00757     vectors[3] = cpl_vector_wrap(nsamples, 
00758             cpl_table_get_data_double((cpl_table*)spc_table, 
00759                 IRPLIB_WLXCORR_COL_OBS)) ;
00760 
00761     /* Scale the signal for a bettre display */
00762     mean1 = cpl_vector_get_mean(vectors[1]) ;
00763     mean3 = cpl_vector_get_mean(vectors[3]) ;
00764     if (fabs(mean3) > 1)
00765         cpl_vector_multiply_scalar(vectors[3], fabs(mean1/mean3)) ;
00766 
00767     cpl_plot_vectors("set grid;set xlabel 'Wavelength (nm)';", title_loc,
00768         "", (const cpl_vector **)vectors, 4);
00769 
00770     /* Unscale the signal */
00771     if (fabs(mean3) > 1)
00772         cpl_vector_multiply_scalar(vectors[3], mean3/mean1) ;
00773 
00774     /* Loop on the brightest lines and zoom on them */
00775     sprintf(title_loc, 
00776 "t '%s - 1-Initial catalog/2-Corrected catalog/3-Observed (ZOOMED)' w lines",
00777         title) ;
00778     title_loc[1023] = (char)0 ;
00779     tmp_vec = cpl_vector_duplicate(vectors[2]) ;
00780     for (i=0 ; i<nblines ; i++) {
00781         /* Find the brightest line */
00782         if ((max = cpl_vector_get_max(tmp_vec)) <= 0.0) break ;
00783         for (j=0 ; i<nsamples ; j++) {
00784             if (cpl_vector_get(tmp_vec, j) == max) break ;
00785         }
00786         if (j-hsize_pix < 0) start_ind = 0 ;
00787         else start_ind = j-hsize_pix ;
00788         if (j+hsize_pix > nsamples-1) stop_ind = nsamples-1 ;
00789         else stop_ind = j+hsize_pix ;
00790         for (j=start_ind ; j<=stop_ind ; j++) cpl_vector_set(tmp_vec, j, 0.0) ;
00791 
00792         sub_vectors = cpl_malloc(4*sizeof(cpl_vector*)) ;
00793         sub_vectors[0]=cpl_vector_extract(vectors[0],start_ind,stop_ind,1);
00794         sub_vectors[1]=cpl_vector_extract(vectors[1],start_ind,stop_ind,1);
00795         sub_vectors[2]=cpl_vector_extract(vectors[2],start_ind,stop_ind,1);
00796         sub_vectors[3]=cpl_vector_extract(vectors[3],start_ind,stop_ind,1);
00797 
00798         cpl_plot_vectors("set grid;set xlabel 'Wavelength (nm)';", title_loc,
00799             "", (const cpl_vector **)sub_vectors, 4);
00800 
00801         cpl_vector_delete(sub_vectors[0]) ;
00802         cpl_vector_delete(sub_vectors[1]) ;
00803         cpl_vector_delete(sub_vectors[2]) ;
00804         cpl_vector_delete(sub_vectors[3]) ;
00805         cpl_free(sub_vectors) ;
00806     }
00807     cpl_vector_delete(tmp_vec) ;
00808     
00809     cpl_vector_unwrap(vectors[0]) ;
00810     cpl_vector_unwrap(vectors[1]) ;
00811     cpl_vector_unwrap(vectors[2]) ;
00812     cpl_vector_unwrap(vectors[3]) ;
00813     cpl_free(vectors) ;
00814 
00815     return 0 ;
00816 }
00817 
00818 /*----------------------------------------------------------------------------*/
00826 /*----------------------------------------------------------------------------*/
00827 int irplib_wlxcorr_catalog_plot(
00828         const cpl_bivector      *   cat,
00829         double                      wmin,
00830         double                      wmax) 
00831 {
00832     int                 start, stop ;
00833     cpl_bivector    *   subcat ;
00834     cpl_vector      *   subcat_x ;
00835     cpl_vector      *   subcat_y ;
00836     const double    *   pwave ;
00837     int                 nvals, nvals_tot ;
00838     int                 i ;
00839 
00840     /* Test entries */
00841     if (cat == NULL) return -1 ;
00842     if (wmax <= wmin) return -1 ;
00843 
00844     /* Initialise */
00845     nvals_tot = cpl_bivector_get_size(cat) ;
00846 
00847     /* Count the nb of values */
00848     pwave = cpl_bivector_get_x_data_const(cat) ;
00849     if (pwave[0] >= wmin) start = 0 ;
00850     else start = -1 ;
00851     if (pwave[nvals_tot-1] <= wmax) stop = nvals_tot-1 ;
00852     else stop = -1 ;
00853     i=0 ;
00854     while ((pwave[i] < wmin) && (i<nvals_tot-1)) i++ ;
00855     start = i ;
00856     i= nvals_tot-1 ;
00857     while ((pwave[i] > wmax) && (i>0)) i-- ;
00858     stop = i ;
00859 
00860     if (start>=stop) {
00861         cpl_msg_error(cpl_func, "Cannot plot the catalog") ;
00862         return -1 ;
00863     }
00864     nvals = start - stop + 1 ;
00865 
00866     /* Create the bivector to plot */
00867     subcat_x = cpl_vector_extract(cpl_bivector_get_x_const(cat),start,stop, 1) ;
00868     subcat_y = cpl_vector_extract(cpl_bivector_get_y_const(cat),start,stop, 1) ;
00869     subcat = cpl_bivector_wrap_vectors(subcat_x, subcat_y) ;
00870 
00871     /* Plot */
00872     if (nvals > 500) {
00873         cpl_plot_bivector(
00874                 "set grid;set xlabel 'Wavelength (nm)';set ylabel 'Emission';",
00875                 "t 'Catalog Spectrum' w lines", "", subcat);
00876     } else {
00877         cpl_plot_bivector(
00878                 "set grid;set xlabel 'Wavelength (nm)';set ylabel 'Emission';",
00879                 "t 'Catalog Spectrum' w impulses", "", subcat);
00880     }
00881     cpl_bivector_unwrap_vectors(subcat) ;
00882     cpl_vector_delete(subcat_x) ;
00883     cpl_vector_delete(subcat_y) ;
00884 
00885     return 0 ;
00886 }
00887    
00890 /*----------------------------------------------------------------------------*/
00905 /*----------------------------------------------------------------------------*/
00906 static void irplib_wlxcorr_estimate(cpl_vector           * vxc,
00907                                     cpl_vector           * model,
00908                                     const cpl_vector     * spectrum,
00909                                     const cpl_bivector   * lines_catalog,
00910                                     const cpl_vector     * conv_kernel,
00911                                     const cpl_polynomial * poly_candi,
00912                                     double                 slitw,
00913                                     double                 fwhm)
00914 {
00915     cpl_errorstate prestate = cpl_errorstate_get();
00916     const int hsize = cpl_vector_get_size(vxc) / 2;
00917 
00918     if (conv_kernel != NULL) {
00919         irplib_wlcalib_fill_spectrum(model, lines_catalog, conv_kernel,
00920                                   poly_candi, hsize);
00921     } else {
00922         const double xtrunc = 0.5 * slitw + 5.0 * fwhm * CPL_MATH_SIG_FWHM;
00923 
00924         irplib_vector_fill_line_spectrum_model(model, NULL, NULL, poly_candi,
00925                                                lines_catalog, slitw, fwhm,
00926                                                xtrunc, 0, CPL_FALSE, CPL_FALSE,
00927                                                NULL);
00928     }
00929 
00930     if (cpl_errorstate_is_equal(prestate))
00931         cpl_vector_correlate(vxc, model, spectrum);
00932 
00933     if (!cpl_errorstate_is_equal(prestate)) {
00934         cpl_vector_fill(vxc, 0.0);
00935 
00936         /* cpl_errorstate_dump(prestate, CPL_FALSE, NULL); */
00937         cpl_errorstate_set(prestate);
00938 
00939     }
00940 
00941     return;
00942 }
00943 
00944 
00945 /*----------------------------------------------------------------------------*/
00955 /*----------------------------------------------------------------------------*/
00956 static cpl_boolean irplib_wlcalib_is_lines(const cpl_vector * wavelengths,
00957                                         const cpl_polynomial * disp1d,
00958                                         int spec_sz,
00959                                         double tol)
00960 {
00961     const int nlines = cpl_vector_get_size(wavelengths);
00962     /* The dispersion on the detector center */
00963     const double dispersion = cpl_polynomial_eval_1d_diff(disp1d,
00964                                                           0.5 * spec_sz + 1.0,
00965                                                           0.5 * spec_sz,
00966                                                           NULL);
00967     const double range = cpl_vector_get(wavelengths, nlines-1)
00968         - cpl_vector_get(wavelengths, 0);
00969 
00970     cpl_ensure(wavelengths != NULL, CPL_ERROR_NULL_INPUT,    CPL_FALSE);
00971     cpl_ensure(disp1d      != NULL, CPL_ERROR_NULL_INPUT,    CPL_FALSE);
00972     cpl_ensure(cpl_polynomial_get_dimension(disp1d) == 1,
00973                CPL_ERROR_ILLEGAL_INPUT, CPL_FALSE);
00974     cpl_ensure(range > 0.0,      CPL_ERROR_ILLEGAL_INPUT, CPL_FALSE);
00975 
00976     return nlines * fabs(dispersion) <= tol * fabs(range) ? CPL_TRUE
00977         : CPL_FALSE;
00978 
00979 }
00980 
00981 /*----------------------------------------------------------------------------*/
00996 /*----------------------------------------------------------------------------*/
00997 static
00998 cpl_error_code irplib_wlcalib_fill_spectrum(cpl_vector           * self,
00999                                          const cpl_bivector   * lines_catalog,
01000                                          const cpl_vector     * conv_kernel,
01001                                          const cpl_polynomial * poly,
01002                                          int                    search_hs)
01003 {
01004 
01005 
01006     const int          size = cpl_vector_get_size(self);
01007     const int          nlines = cpl_bivector_get_size(lines_catalog);
01008     const cpl_vector * xlines  = cpl_bivector_get_x_const(lines_catalog);
01009     const double     * dxlines = cpl_vector_get_data_const(xlines);
01010     cpl_bivector     * sub_cat ;
01011     cpl_vector       * sub_cat_x;
01012     cpl_vector       * sub_cat_y;
01013     cpl_vector       * wl_limits;
01014     double             wave_min, wave_max;
01015     int                wave_min_id, wave_max_id;
01016     int                nsub;
01017     int                error;
01018 
01019     cpl_ensure_code(self          != NULL, CPL_ERROR_NULL_INPUT);
01020     cpl_ensure_code(lines_catalog != NULL, CPL_ERROR_NULL_INPUT);
01021     cpl_ensure_code(conv_kernel   != NULL, CPL_ERROR_NULL_INPUT);
01022     cpl_ensure_code(poly          != NULL, CPL_ERROR_NULL_INPUT);
01023     cpl_ensure_code(size > 0,              CPL_ERROR_ILLEGAL_INPUT);
01024 
01025 
01026     /* Resample the spectrum */
01027     wl_limits = cpl_vector_new(size + 1);
01028     cpl_vector_fill_polynomial(wl_limits, poly, 0.5 - search_hs, 1);
01029 
01030     /* The spectrum wavelength bounds */
01031     wave_min = cpl_vector_get(wl_limits, 0);
01032     wave_max = cpl_vector_get(wl_limits, size);
01033 
01034     /* Find the 1st line */
01035     wave_min_id = cpl_vector_find(xlines, wave_min);
01036     /* The first line must be less than or equal to wave_min */
01037     if (dxlines[wave_min_id] > wave_min) wave_min_id--;
01038 
01039     if (wave_min_id < 0) {
01040         cpl_vector_delete(wl_limits);
01041         return cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
01042                                            __FILE__, __LINE__, "The %d-line "
01043                                            "catalogue only has lines above %g",
01044                                            nlines, wave_min);
01045     }
01046 
01047     /* Find the last line */
01048     wave_max_id = cpl_vector_find(xlines, wave_max);
01049     /* The last line must be greater than or equal to wave_max */
01050     if (dxlines[wave_max_id] < wave_max) wave_max_id++;
01051 
01052     if (wave_max_id == nlines) {
01053         cpl_vector_delete(wl_limits);
01054         return cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
01055                                            __FILE__, __LINE__, "The %d-line "
01056                                            "catalogue only has lines below %g",
01057                                            nlines, wave_max);
01058     }
01059 
01060     /* Checking the wavelength range at this point via the indices also
01061        verifies that they were not found using non-increasing wavelengths */
01062     nsub = 1 + wave_max_id - wave_min_id;
01063     cpl_ensure_code(nsub > 1, CPL_ERROR_ILLEGAL_INPUT);
01064 
01065     /* Wrap a new bivector around the relevant part of the catalog */
01066     /* The data is _not_ modified */
01067     sub_cat_x = cpl_vector_wrap(nsub, wave_min_id + (double*)dxlines);
01068     sub_cat_y = cpl_vector_wrap(nsub, wave_min_id + (double*)
01069                                 cpl_bivector_get_y_data_const(lines_catalog));
01070     sub_cat = cpl_bivector_wrap_vectors(sub_cat_x, sub_cat_y);
01071 
01072     /* High resolution catalog */
01073     error = irplib_wlxcorr_signal_resample(self, wl_limits, sub_cat);
01074 
01075     cpl_vector_delete(wl_limits);
01076     cpl_bivector_unwrap_vectors(sub_cat);
01077     (void)cpl_vector_unwrap(sub_cat_x);
01078     (void)cpl_vector_unwrap(sub_cat_y);
01079 
01080     cpl_ensure_code(!error, CPL_ERROR_ILLEGAL_INPUT);
01081 
01082     /* Smooth the instrument resolution */
01083     cpl_ensure_code(!irplib_wlxcorr_convolve(self, conv_kernel),
01084                     cpl_error_get_code());
01085 
01086     return CPL_ERROR_NONE;
01087 }
01088 
01089 
01090 /*----------------------------------------------------------------------------*/
01100 /*----------------------------------------------------------------------------*/
01101 static int irplib_wlxcorr_signal_resample(
01102         cpl_vector          *   resampled, 
01103         const cpl_vector    *   xbounds,
01104         const cpl_bivector  *   hires)
01105 {
01106     const int           hrsize = cpl_bivector_get_size(hires);
01107     const cpl_vector*   xhires ;
01108     const cpl_vector*   yhires ;
01109     const double    *   pxhires ;
01110     const double    *   pyhires ;
01111     const double    *   pxbounds ;
01112     cpl_vector      *   ybounds ;
01113     cpl_bivector    *   boundary ;
01114     double          *   pybounds ;
01115     double          *   presampled ;
01116     int                 nsamples ;
01117     int                 i, itt ;
01118    
01119     /* Test entries */
01120     if ((!resampled) || (!xbounds) || (!hires)) return -1 ;
01121 
01122     /* Initialise */
01123     nsamples = cpl_vector_get_size(resampled) ;
01124 
01125     /* Initialise */
01126     presampled = cpl_vector_get_data(resampled) ;
01127     pxbounds = cpl_vector_get_data_const(xbounds) ;
01128     xhires = cpl_bivector_get_x_const(hires) ;
01129     yhires = cpl_bivector_get_y_const(hires) ;
01130     pxhires = cpl_vector_get_data_const(xhires) ;
01131     pyhires = cpl_vector_get_data_const(yhires) ;
01132     
01133     /* Create a new vector */
01134     ybounds = cpl_vector_new(cpl_vector_get_size(xbounds)) ;
01135     boundary = cpl_bivector_wrap_vectors((cpl_vector*)xbounds,ybounds) ;
01136     pybounds = cpl_vector_get_data(ybounds) ;
01137 
01138     /* Test entries */
01139     if (cpl_bivector_get_size(boundary) != nsamples + 1) {
01140         cpl_bivector_unwrap_vectors(boundary) ;
01141         cpl_vector_delete(ybounds) ;
01142         return -1 ;
01143     }
01144 
01145     /* Get the ind  */
01146     itt = cpl_vector_find(xhires, pxbounds[0]);
01147 
01148     /* Interpolate the signal */
01149     if (cpl_bivector_interpolate_linear(boundary, hires)) {
01150         cpl_bivector_unwrap_vectors(boundary) ;
01151         cpl_vector_delete(ybounds) ;
01152         return -1 ;
01153     }
01154 
01155     /* At this point itt most likely points to element just below
01156        pxbounds[0] */
01157     while (pxhires[itt] < pxbounds[0]) itt++;
01158 
01159     for (i=0; i < nsamples; i++) {
01160         /* The i'th signal is the weighted average of the two interpolated
01161            signals at the pixel boundaries and those table signals in
01162            between */
01163 
01164         double xlow  = pxbounds[i];
01165         double x     = pxhires[itt];
01166 
01167         if (x > pxbounds[i+1]) x = pxbounds[i+1];
01168         /* Contribution from interpolated value at wavelength at lower pixel
01169            boundary */
01170         presampled[i] = pybounds[i] * (x - xlow);
01171 
01172         /* Contribution from table values in between pixel boundaries */
01173         while ((pxhires[itt] < pxbounds[i+1]) && (itt < hrsize)) {
01174             const double xprev = x;
01175             x = pxhires[itt+1];
01176             if (x > pxbounds[i+1]) x = pxbounds[i+1];
01177             presampled[i] += pyhires[itt] * (x - xlow);
01178             xlow = xprev;
01179             itt++;
01180         }
01181 
01182         /* Contribution from interpolated value at wavelength at upper pixel
01183            boundary */
01184         presampled[i] += pybounds[i+1] * (pxbounds[i+1] - xlow);
01185 
01186         /* Compute average by dividing integral by length of pixel range
01187            (the factor 2 comes from the contributions) */
01188         presampled[i] /= 2 * (pxbounds[i+1] - pxbounds[i]);
01189     }
01190     cpl_bivector_unwrap_vectors(boundary) ;
01191     cpl_vector_delete(ybounds) ;
01192     return 0 ;
01193 }
01194 
01195 
01196 
01197 /*----------------------------------------------------------------------------*/
01218 /*----------------------------------------------------------------------------*/
01219 static cpl_error_code cpl_vector_fill_lss_profile_symmetric(cpl_vector * self,
01220                                                             double  slitw,
01221                                                             double  fwhm)
01222 {
01223 
01224     const double sigma = fwhm * CPL_MATH_SIG_FWHM;
01225     const int    n     = cpl_vector_get_size(self);
01226     int          i;
01227 
01228 
01229     cpl_ensure_code(self != NULL, CPL_ERROR_NULL_INPUT);
01230     cpl_ensure_code(slitw > 0.0,  CPL_ERROR_ILLEGAL_INPUT);
01231     cpl_ensure_code(fwhm  > 0.0,  CPL_ERROR_ILLEGAL_INPUT);
01232 
01233     /* Cannot fail now */
01234 
01235     /* Special case for i = 0 */
01236     (void)cpl_vector_set(self, 0,
01237                          (irplib_erf_antideriv(0.5*slitw + 0.5, sigma) -
01238                           irplib_erf_antideriv(0.5*slitw - 0.5, sigma)) / slitw);
01239 
01240     for (i = 1; i < n; i++) {
01241         /* FIXME: Reuse two irplib_erf_antideriv() calls from previous value */
01242         const double x1p = i + 0.5*slitw + 0.5;
01243         const double x1n = i - 0.5*slitw + 0.5;
01244         const double x0p = i + 0.5*slitw - 0.5;
01245         const double x0n = i - 0.5*slitw - 0.5;
01246         const double val = 0.5/slitw *
01247             (irplib_erf_antideriv(x1p, sigma) - irplib_erf_antideriv(x1n, sigma) -
01248              irplib_erf_antideriv(x0p, sigma) + irplib_erf_antideriv(x0n, sigma));
01249         (void)cpl_vector_set(self, i, val);
01250     }
01251 
01252     return CPL_ERROR_NONE;
01253 }

Generated on 1 Mar 2011 for DETMON Pipeline Reference Manual by  doxygen 1.6.1