irplib_wlxcorr.c

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

Generated on 7 Mar 2012 for DETMON Pipeline Reference Manual by  doxygen 1.6.1