irplib_wavecal.c

00001 /* $Id: irplib_wavecal.c,v 1.46 2012/03/02 09:01:04 amodigli Exp $
00002  *
00003  * This file is part of the IRPLIB Pipeline
00004  * Copyright (C) 2002,2003 European Southern Observatory
00005  *
00006  * This program is free software; you can redistribute it and/or modify
00007  * it under the terms of the GNU General Public License as published by
00008  * the Free Software Foundation; either version 2 of the License, or
00009  * (at your option) any later version.
00010  *
00011  * This program is distributed in the hope that it will be useful,
00012  * but WITHOUT ANY WARRANTY; without even the implied warranty of
00013  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00014  * GNU General Public License for more details.
00015  *
00016  * You should have received a copy of the GNU General Public License
00017  * along with this program; if not, write to the Free Software
00018  * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02111-1307  USA
00019  */
00020 
00021 /*
00022  * $Author: amodigli $
00023  * $Date: 2012/03/02 09:01:04 $
00024  * $Revision: 1.46 $
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 /* TEMPORARY SUPPORT OF CPL 5.x */
00037 #include <cpl.h>
00038 
00039 #ifndef CPL_SIZE_FORMAT
00040 #define CPL_SIZE_FORMAT "d"
00041 #define cpl_size int
00042 #endif
00043 /* END TEMPORARY SUPPORT OF CPL 5.x */
00044 
00045 #include "irplib_wavecal_impl.h"
00046 
00047 /* Needed for irplib_errorstate_dump_debug() */
00048 #include "irplib_utils.h"
00049 
00050 #include <string.h>
00051 #include <math.h>
00052 
00053 #ifdef HAVE_GSL
00054 #include <gsl/gsl_multimin.h>
00055 #endif
00056 
00057 /*-----------------------------------------------------------------------------
00058                                Private types
00059  -----------------------------------------------------------------------------*/
00060 
00061 typedef struct {
00062 
00063     const cpl_vector * observed;
00064     cpl_polynomial   * disp1d;
00065     cpl_vector       * spectrum;
00066     irplib_base_spectrum_model * param;
00067     cpl_error_code  (* filler)(cpl_vector *, const cpl_polynomial *,
00068                                irplib_base_spectrum_model *, int);
00069     cpl_vector       * vxc;
00070     double             xc;
00071     int                maxxc;
00072     double             mxc;
00073     cpl_polynomial   * mdisp;
00074     int                ishift;
00075 
00076 } irplib_multimin;
00077 
00078 /*-----------------------------------------------------------------------------
00079                                Defines
00080  -----------------------------------------------------------------------------*/
00081 
00082 #ifndef inline
00083 #define inline /* inline */
00084 #endif
00085 
00086 #define IRPLIB_MAX(A,B) ((A) > (B) ? (A) : (B))
00087 #define IRPLIB_MIN(A,B) ((A) < (B) ? (A) : (B))
00088 
00089 /*-----------------------------------------------------------------------------
00090                                    Private functions
00091  -----------------------------------------------------------------------------*/
00092 
00093 #ifdef HAVE_GSL
00094 static double irplib_gsl_correlation(const gsl_vector *, void *);
00095 #endif
00096 
00097 /*----------------------------------------------------------------------------*/
00101 /*----------------------------------------------------------------------------*/
00102 
00106 /*----------------------------------------------------------------------------*/
00114 /*----------------------------------------------------------------------------*/
00115 int irplib_bivector_count_positive(const cpl_bivector * self,
00116                                   double               x_min,
00117                                   double               x_max)
00118 {
00119 
00120     const int      nself = cpl_bivector_get_size(self);
00121     const double * px    = cpl_bivector_get_x_data_const(self);
00122     const double * py    = cpl_bivector_get_y_data_const(self);
00123     int            npos  = 0;
00124     int            i     = 0;
00125 
00126     cpl_ensure(self != NULL, CPL_ERROR_NULL_INPUT, -1);
00127     cpl_ensure(x_min <= x_max, CPL_ERROR_ILLEGAL_INPUT, -2);
00128     
00129     /* FIXME: Use cpl_vector_find() */
00130     while (i < nself && px[i] < x_min) i++;
00131     while (i < nself && px[i] < x_max)
00132         if (py[i++] > 0) npos++;
00133 
00134     return npos;
00135 }
00136 
00137 /*----------------------------------------------------------------------------*/
00147 /*----------------------------------------------------------------------------*/
00148 cpl_error_code irplib_polynomial_fit_2d_dispersion(cpl_polynomial * self,
00149                                                    const cpl_image * imgwave,
00150                                                    int fitdeg, double * presid)
00151 {
00152 
00153     const int        nx = cpl_image_get_size_x(imgwave);
00154     const int        ny = cpl_image_get_size_y(imgwave);
00155     const int        nbad = cpl_image_count_rejected(imgwave);
00156     const int        nsamp = nx * ny - nbad;
00157     cpl_matrix     * xy_pos;
00158     double         * xdata;
00159     double         * ydata;
00160     cpl_vector     * wlen;
00161     double         * dwlen;
00162     const cpl_size   nfitdeg = (cpl_size)fitdeg;
00163     int i, j;
00164     int k = 0;
00165 
00166     cpl_ensure_code(self    != NULL, CPL_ERROR_NULL_INPUT);
00167     cpl_ensure_code(imgwave != NULL, CPL_ERROR_NULL_INPUT);
00168     cpl_ensure_code(presid  != NULL, CPL_ERROR_NULL_INPUT);
00169     cpl_ensure_code(fitdeg > 0,      CPL_ERROR_ILLEGAL_INPUT);
00170 
00171     cpl_ensure_code(cpl_polynomial_get_dimension(self) == 2,
00172                     CPL_ERROR_ILLEGAL_INPUT);
00173 
00174     xy_pos = cpl_matrix_new(2, nsamp);
00175     xdata = cpl_matrix_get_data(xy_pos);
00176     ydata = xdata + nsamp;
00177 
00178     dwlen = (double*)cpl_malloc(nsamp * sizeof(double));
00179     wlen = cpl_vector_wrap(nsamp, dwlen);
00180 
00181     for (i=1; i <= nx; i++) {
00182         for (j=1; j <= ny; j++) {
00183             int is_bad;
00184             const double value = cpl_image_get(imgwave, i, j, &is_bad);
00185             if (!is_bad) {
00186                 xdata[k] = i;
00187                 ydata[k] = j;
00188                 dwlen[k] = value;
00189                 k++;
00190             }
00191         }
00192     }
00193 
00194     cpl_msg_info(cpl_func, "Fitting 2D polynomial to %d X %d image, ignoring "
00195                  "%d poorly calibrated pixels", nx, ny, nbad);
00196 
00197     if (cpl_polynomial_fit(self, xy_pos, NULL, wlen, NULL, CPL_FALSE, NULL,
00198                            &nfitdeg) == CPL_ERROR_NONE && presid != NULL) {
00199         cpl_vector_fill_polynomial_fit_residual(wlen, wlen, NULL, self, xy_pos,
00200                                                 NULL);
00201         *presid = cpl_vector_product(wlen, wlen)/nsamp;
00202     }
00203     cpl_matrix_delete(xy_pos);
00204     cpl_vector_delete(wlen);
00205 
00206     cpl_ensure_code(k == nsamp, CPL_ERROR_UNSPECIFIED);
00207 
00208     return CPL_ERROR_NONE;
00209 }
00210 
00211 
00212 /*----------------------------------------------------------------------------*/
00230 /*----------------------------------------------------------------------------*/
00231 cpl_error_code
00232 irplib_polynomial_find_1d_from_correlation(cpl_polynomial * self,
00233                                            int maxdeg,
00234                                            const cpl_vector * obs,
00235                                            irplib_base_spectrum_model * model,
00236                                            cpl_error_code (* filler)
00237                                            (cpl_vector *,
00238                                             const cpl_polynomial *,
00239                                             irplib_base_spectrum_model *, int),
00240                                            double pixtol,
00241                                            double pixstep,
00242                                            int hsize,
00243                                            int maxite,
00244                                            double * pxc)
00245 {
00246 
00247 #ifdef HAVE_GSL
00248     const gsl_multimin_fminimizer_type * T = gsl_multimin_fminimizer_nmsimplex;
00249     gsl_multimin_fminimizer * minimizer;
00250     gsl_multimin_function my_func;
00251     irplib_multimin data;
00252     gsl_vector * dispgsl;
00253     gsl_vector * stepsize;
00254     gsl_vector * dispprev;
00255     int status = GSL_CONTINUE;
00256     const int nobs = cpl_vector_get_size(obs);
00257     const cpl_size nfit = maxdeg + 1;
00258     cpl_errorstate prestate = cpl_errorstate_get();
00259     /* Convert pixel step to wavelength step on detector center */
00260     const double wlstep =
00261         cpl_polynomial_eval_1d_diff(self, 0.5 * (nobs + pixstep),
00262                                     0.5 * (nobs - pixstep), NULL);
00263     double wlstepi = wlstep;
00264     double size;
00265     int iter;
00266     cpl_size i;
00267 
00268 #endif
00269 
00270     cpl_ensure_code(self   != NULL, CPL_ERROR_NULL_INPUT);
00271     cpl_ensure_code(obs    != NULL, CPL_ERROR_NULL_INPUT);
00272     cpl_ensure_code(model  != NULL, CPL_ERROR_NULL_INPUT);
00273     cpl_ensure_code(filler != NULL, CPL_ERROR_NULL_INPUT);
00274     cpl_ensure_code(pxc    != NULL, CPL_ERROR_NULL_INPUT);
00275 
00276     cpl_ensure_code(cpl_polynomial_get_dimension(self) == 1,
00277                     CPL_ERROR_ILLEGAL_INPUT);
00278 
00279     cpl_ensure_code(cpl_polynomial_get_degree(self) > 0,
00280                     CPL_ERROR_ILLEGAL_INPUT);
00281 
00282     cpl_ensure_code(maxdeg  >=  0, CPL_ERROR_ILLEGAL_INPUT);
00283     cpl_ensure_code(pixtol  > 0.0, CPL_ERROR_ILLEGAL_INPUT);
00284     cpl_ensure_code(pixstep > 0.0, CPL_ERROR_ILLEGAL_INPUT);
00285     cpl_ensure_code(hsize   >=  0, CPL_ERROR_ILLEGAL_INPUT);
00286     cpl_ensure_code(maxite  >=  0, CPL_ERROR_ILLEGAL_INPUT);
00287 
00288 #ifndef HAVE_GSL
00289     return cpl_error_set_message(cpl_func, CPL_ERROR_UNSUPPORTED_MODE,
00290                                  "GSL is not available");
00291 #else
00292 
00293     minimizer = gsl_multimin_fminimizer_alloc(T, (size_t)nfit);
00294 
00295     cpl_ensure_code(minimizer != NULL, CPL_ERROR_ILLEGAL_OUTPUT);
00296        
00297     dispgsl  = gsl_vector_alloc((size_t)nfit);
00298     stepsize = gsl_vector_alloc((size_t)nfit);
00299     dispprev = gsl_vector_alloc((size_t)nfit);
00300 
00301     for (i=0; i < nfit; i++) {
00302         const double value = cpl_polynomial_get_coeff(self, &i);
00303         gsl_vector_set(dispgsl, (size_t)i, value);
00304         gsl_vector_set(stepsize, (size_t)i, wlstepi);
00305         wlstepi /= (double)nobs;
00306     }
00307 
00308     my_func.n = nfit;
00309     my_func.f = &irplib_gsl_correlation;
00310     my_func.params = (void *)(&data);
00311 
00312     data.observed = obs;
00313     data.disp1d   = self;
00314     data.spectrum = cpl_vector_new(nobs + 2 * hsize);
00315     data.vxc      = cpl_vector_new(1 + 2 * hsize);
00316     data.param    = model;
00317     data.filler   = filler;
00318     data.maxxc    = 0; /* Output */
00319     data.ishift   = 0; /* Output */
00320     data.mxc      = -1.0; /* Output */
00321     data.mdisp    = NULL; /* Output */
00322 
00323     gsl_multimin_fminimizer_set (minimizer, &my_func, dispgsl, stepsize);
00324 
00325     for (iter = 0; status == GSL_CONTINUE && iter < maxite; iter++) {
00326 
00327         const double fprev = minimizer->fval;
00328 
00329         gsl_vector_memcpy(dispprev, minimizer->x);
00330         status = gsl_multimin_fminimizer_iterate(minimizer);
00331 
00332         if (status || !cpl_errorstate_is_equal(prestate)) break;
00333 
00334         size = gsl_multimin_fminimizer_size (minimizer);
00335         status = gsl_multimin_test_size (size, pixtol);
00336      
00337         if (status == GSL_SUCCESS) {
00338             cpl_msg_debug(cpl_func, "converged to minimum at");
00339 
00340             if (nfit == 0) {
00341                 cpl_msg_debug(cpl_func, "%5d %g df() = %g size = %g", 
00342                               iter,
00343                               gsl_vector_get (minimizer->x, 0)
00344                               - gsl_vector_get (dispprev, 0), 
00345                               minimizer->fval - fprev, size);
00346             } else if (nfit == 1) {
00347                 cpl_msg_debug(cpl_func, "%5d %g %g df() = %g size = %g", 
00348                               iter,
00349                               gsl_vector_get (minimizer->x, 0)
00350                               - gsl_vector_get (dispprev, 0), 
00351                               gsl_vector_get (minimizer->x, 1)
00352                               - gsl_vector_get (dispprev, 1), 
00353                               minimizer->fval - fprev, size);
00354             } else {
00355                 cpl_msg_debug(cpl_func, "%5d %g %g %g df() = %g size = %g", 
00356                               iter,
00357                               gsl_vector_get (minimizer->x, 0)
00358                               - gsl_vector_get (dispprev, 0), 
00359                               gsl_vector_get (minimizer->x, 1)
00360                               - gsl_vector_get (dispprev, 1), 
00361                               gsl_vector_get (minimizer->x, 2)
00362                               - gsl_vector_get (dispprev, 2), 
00363                               minimizer->fval - fprev, size);
00364             }
00365         }
00366     }
00367 
00368     if (status == GSL_SUCCESS && cpl_errorstate_is_equal(prestate)) {
00369         if (data.mxc > -minimizer->fval) {
00370             *pxc = data.mxc;
00371             cpl_msg_warning(cpl_func, "Local maximum: %g(%d) > %g",
00372                             data.mxc, data.ishift, -minimizer->fval);
00373             cpl_polynomial_shift_1d(data.mdisp, 0, (double)data.ishift);
00374             cpl_polynomial_copy(self, data.mdisp);
00375             status = GSL_CONTINUE;
00376         } else {
00377             *pxc = -minimizer->fval;
00378             for (i=0; i < nfit; i++) {
00379                 const double value = gsl_vector_get(minimizer->x, i);
00380                 cpl_polynomial_set_coeff(self, &i, value);
00381             }
00382         }
00383     }
00384 
00385     cpl_vector_delete(data.spectrum);
00386     cpl_vector_delete(data.vxc);
00387     cpl_polynomial_delete(data.mdisp);
00388     gsl_multimin_fminimizer_free(minimizer);
00389     gsl_vector_free(dispgsl);
00390     gsl_vector_free(dispprev);
00391     gsl_vector_free(stepsize);
00392 
00393     cpl_ensure_code(status != GSL_CONTINUE, CPL_ERROR_CONTINUE);
00394     cpl_ensure_code(status == GSL_SUCCESS, CPL_ERROR_DATA_NOT_FOUND);
00395     cpl_ensure_code(cpl_errorstate_is_equal(prestate), cpl_error_get_code());
00396 
00397     return CPL_ERROR_NONE;
00398 #endif
00399 }
00400 
00401 
00402 /*----------------------------------------------------------------------------*/
00431 /*----------------------------------------------------------------------------*/
00432 cpl_error_code
00433 irplib_vector_fill_line_spectrum(cpl_vector * self,
00434                                  const cpl_polynomial * disp,
00435                                  irplib_base_spectrum_model * lsslamp,
00436                                  int hsize)
00437 {
00438 
00439     irplib_line_spectrum_model * arclamp
00440         = (irplib_line_spectrum_model *)lsslamp;
00441     cpl_error_code error;
00442 
00443     cpl_ensure_code(arclamp != NULL, CPL_ERROR_NULL_INPUT);
00444 
00445     arclamp->cost++;
00446 
00447     error = irplib_vector_fill_line_spectrum_model(self,
00448                                                    arclamp->linepix,
00449                                                    arclamp->erftmp,
00450                                                    disp,
00451                                                    arclamp->lines,
00452                                                    arclamp->wslit,
00453                                                    arclamp->wfwhm,
00454                                                    arclamp->xtrunc,
00455                                                    hsize, CPL_FALSE, CPL_FALSE,
00456                                                    &(arclamp->ulines));
00457     cpl_ensure_code(!error, error);
00458 
00459     arclamp->xcost++;
00460 
00461     return CPL_ERROR_NONE;
00462 }
00463 
00464 /*----------------------------------------------------------------------------*/
00478 /*----------------------------------------------------------------------------*/
00479 cpl_error_code
00480 irplib_vector_fill_logline_spectrum(cpl_vector * self,
00481                                     const cpl_polynomial * disp,
00482                                     irplib_base_spectrum_model * lsslamp,
00483                                     int hsize)
00484 {
00485 
00486     irplib_line_spectrum_model * arclamp
00487         = (irplib_line_spectrum_model *)lsslamp;
00488     cpl_error_code error;
00489 
00490     cpl_ensure_code(arclamp != NULL, CPL_ERROR_NULL_INPUT);
00491 
00492     arclamp->cost++;
00493 
00494     error = irplib_vector_fill_line_spectrum_model(self,
00495                                                    arclamp->linepix,
00496                                                    arclamp->erftmp,
00497                                                    disp,
00498                                                    arclamp->lines,
00499                                                    arclamp->wslit,
00500                                                    arclamp->wfwhm,
00501                                                    arclamp->xtrunc,
00502                                                    hsize, CPL_FALSE, CPL_TRUE,
00503                                                    &(arclamp->ulines));
00504     cpl_ensure_code(!error, error);
00505 
00506     arclamp->xcost++;
00507 
00508     return CPL_ERROR_NONE;
00509 }
00510 
00511 
00512 /*----------------------------------------------------------------------------*/
00526 /*----------------------------------------------------------------------------*/
00527 cpl_error_code
00528 irplib_vector_fill_line_spectrum_fast(cpl_vector * self,
00529                                       const cpl_polynomial * disp,
00530                                       irplib_base_spectrum_model * lsslamp,
00531                                       int hsize)
00532 {
00533 
00534     irplib_line_spectrum_model * arclamp
00535         = (irplib_line_spectrum_model *)lsslamp;
00536     cpl_error_code error;
00537 
00538     cpl_ensure_code(arclamp != NULL, CPL_ERROR_NULL_INPUT);
00539 
00540     arclamp->cost++;
00541 
00542     error = irplib_vector_fill_line_spectrum_model(self,
00543                                                    arclamp->linepix,
00544                                                    arclamp->erftmp,
00545                                                    disp,
00546                                                    arclamp->lines,
00547                                                    arclamp->wslit,
00548                                                    arclamp->wfwhm,
00549                                                    arclamp->xtrunc,
00550                                                    hsize, CPL_TRUE, CPL_FALSE,
00551                                                    &(arclamp->ulines));
00552     cpl_ensure_code(!error, error);
00553 
00554     arclamp->xcost++;
00555 
00556     return CPL_ERROR_NONE;
00557 }
00558 
00559 /*----------------------------------------------------------------------------*/
00573 /*----------------------------------------------------------------------------*/
00574 cpl_error_code
00575 irplib_vector_fill_logline_spectrum_fast(cpl_vector * self,
00576                                          const cpl_polynomial * disp,
00577                                          irplib_base_spectrum_model * lsslamp,
00578                                          int hsize)
00579 {
00580 
00581     irplib_line_spectrum_model * arclamp
00582         = (irplib_line_spectrum_model *)lsslamp;
00583     cpl_error_code error;
00584 
00585     cpl_ensure_code(arclamp != NULL, CPL_ERROR_NULL_INPUT);
00586 
00587     arclamp->cost++;
00588 
00589     error = irplib_vector_fill_line_spectrum_model(self,
00590                                                    arclamp->linepix,
00591                                                    arclamp->erftmp,
00592                                                    disp,
00593                                                    arclamp->lines,
00594                                                    arclamp->wslit,
00595                                                    arclamp->wfwhm,
00596                                                    arclamp->xtrunc,
00597                                                    hsize, CPL_TRUE, CPL_TRUE,
00598                                                    &(arclamp->ulines));
00599     cpl_ensure_code(!error, error);
00600 
00601     arclamp->xcost++;
00602 
00603     return CPL_ERROR_NONE;
00604 }
00605 
00606 /*----------------------------------------------------------------------------*/
00617 /*----------------------------------------------------------------------------*/
00618 cpl_error_code irplib_plot_spectrum_and_model(const cpl_vector * self,
00619                                               const cpl_polynomial * disp1d,
00620                                               irplib_base_spectrum_model * model,
00621                                               cpl_error_code (* filler)
00622                                               (cpl_vector *,
00623                                                const cpl_polynomial *,
00624                                                irplib_base_spectrum_model *,
00625                                                int))
00626 {
00627 
00628     cpl_errorstate prestate = cpl_errorstate_get();
00629     cpl_vector * wl;
00630     cpl_vector * spectrum;
00631     cpl_vector * vxc;
00632     const int len = cpl_vector_get_size(self);
00633     double maxval, xc;
00634     int ixc;
00635     int error = 0;
00636 
00637     cpl_ensure_code(self   != NULL, CPL_ERROR_NULL_INPUT);
00638     cpl_ensure_code(disp1d != NULL, CPL_ERROR_NULL_INPUT);
00639     cpl_ensure_code(model  != NULL, CPL_ERROR_NULL_INPUT);
00640     cpl_ensure_code(filler != NULL, CPL_ERROR_NULL_INPUT);
00641 
00642     cpl_ensure_code(cpl_polynomial_get_dimension(disp1d) == 1,
00643                     CPL_ERROR_ILLEGAL_INPUT);
00644 
00645     cpl_ensure_code(cpl_polynomial_get_degree(disp1d) > 0,
00646                     CPL_ERROR_ILLEGAL_INPUT);
00647 
00648     wl = cpl_vector_new(len);
00649     spectrum = cpl_vector_new(len);
00650     vxc = cpl_vector_new(1);
00651 
00652     error |= (int)cpl_vector_fill_polynomial(wl, disp1d, 1.0, 1.0);
00653     error |= filler(spectrum, disp1d, model, 0);
00654 
00655     ixc = cpl_vector_correlate(vxc, self, spectrum);
00656     xc = cpl_vector_get(vxc, ixc);
00657 
00658     maxval = cpl_vector_get_max(spectrum);
00659     if (maxval != 0.0) 
00660         error |= cpl_vector_multiply_scalar(spectrum,
00661                                              cpl_vector_get_max(self)/maxval);
00662     if (!error) {
00663         const cpl_vector * spair[] = {wl, self, spectrum};
00664         char * pre = cpl_sprintf("set grid;set xlabel 'Wavelength (%g -> %g)'; "
00665                                  "set ylabel 'Intensity';", cpl_vector_get(wl, 0),
00666                                  cpl_vector_get(wl, len-1));
00667         char * title = cpl_sprintf("t 'Observed and modelled spectra (%d pixel "
00668                                    "XC=%g) ' w linespoints", len, xc);
00669 
00670         (void)cpl_plot_vectors(pre, title, "", spair, 3);
00671         cpl_free(pre);
00672         cpl_free(title);
00673     }
00674 
00675     cpl_vector_delete(wl);
00676     cpl_vector_delete(spectrum);
00677     cpl_vector_delete(vxc);
00678 
00679     cpl_errorstate_set(prestate);
00680 
00681     return CPL_ERROR_NONE;
00682 }
00683 
00684 /*----------------------------------------------------------------------------*/
00702 /*----------------------------------------------------------------------------*/
00703 cpl_error_code
00704 irplib_bivector_find_shift_from_correlation(cpl_bivector * self,
00705                                             const cpl_polynomial * disp,
00706                                             const cpl_vector * obs,
00707                                             irplib_base_spectrum_model * model,
00708                                             cpl_error_code (*filler)
00709                                             (cpl_vector *,
00710                                              const cpl_polynomial *,
00711                                              irplib_base_spectrum_model *, int),
00712                                             int hsize,
00713                                             cpl_boolean doplot,
00714                                             double *pxc)
00715 {
00716 
00717     const int      nobs   = cpl_vector_get_size(obs);
00718     const int      nmodel = 2 * hsize + nobs;
00719     cpl_vector   * xself = cpl_bivector_get_x(self);
00720     cpl_vector   * yself = cpl_bivector_get_y(self);
00721     cpl_vector   * mspec1d;
00722     cpl_vector   * xcorr;
00723     cpl_error_code error = CPL_ERROR_NONE;
00724     double         xcprev, xcnext;
00725     int            ixc, imax = 0;
00726     int i;
00727 
00728     cpl_ensure_code(self   != NULL, CPL_ERROR_NULL_INPUT);
00729     cpl_ensure_code(disp   != NULL, CPL_ERROR_NULL_INPUT);
00730     cpl_ensure_code(obs    != NULL, CPL_ERROR_NULL_INPUT);
00731     cpl_ensure_code(model  != NULL, CPL_ERROR_NULL_INPUT);
00732     cpl_ensure_code(filler != NULL, CPL_ERROR_NULL_INPUT);
00733     cpl_ensure_code(hsize  >  0,    CPL_ERROR_ILLEGAL_INPUT);
00734 
00735     mspec1d = cpl_vector_new(nmodel);
00736 
00737     if (filler(mspec1d, disp, model, hsize)) {
00738         cpl_vector_delete(mspec1d);
00739         return cpl_error_set_where(cpl_func);
00740     }
00741 
00742     /* Should not be able to fail now */
00743     xcorr = cpl_vector_new(1 + 2 * hsize);
00744     ixc = cpl_vector_correlate(xcorr, mspec1d, obs);
00745 
00746 #ifdef IRPLIB_SPC_DUMP
00747     /* Need irplib_wavecal.c rev. 1.12 through 1.15 */
00748     irplib_polynomial_dump_corr_step(disp, xcorr, "Shift");
00749 #endif
00750 
00751     cpl_vector_delete(mspec1d);
00752 
00753     /* Find local maxima. */
00754     /* FIXME(?): Also include stationary points */
00755     i = 0;
00756     xcprev = cpl_vector_get(xcorr, i);
00757     xcnext = cpl_vector_get(xcorr, i+1);
00758 
00759     if (xcprev >= xcnext) {
00760         /* 1st data point is an extreme */
00761         /* FIXME: This could also be an error, recoverable by caller by
00762            increasing hsize */
00763         imax++;
00764 
00765         cpl_vector_set(xself, 0, i - hsize);
00766         cpl_vector_set(yself, 0, xcprev);
00767 
00768     }
00769 
00770     for (i = 1; i < 2 * hsize; i++) {
00771         const double xc = xcnext;
00772         xcnext = cpl_vector_get(xcorr, i+1);
00773         if (xc >= xcprev && xc >= xcnext) {
00774             /* Found (local) maximum at shift i - hsize */
00775             int j;
00776 
00777             imax++;
00778 
00779             if (cpl_bivector_get_size(self) < imax) {
00780                 cpl_vector_set_size(xself, imax);
00781                 cpl_vector_set_size(yself, imax);
00782             }
00783 
00784             for (j = imax-1; j > 0; j--) {
00785                 if (xc <= cpl_vector_get(yself, j-1)) break;
00786                 cpl_vector_set(xself, j, cpl_vector_get(xself, j-1));
00787                 cpl_vector_set(yself, j, cpl_vector_get(yself, j-1));
00788             }
00789             cpl_vector_set(xself, j, i - hsize);
00790             cpl_vector_set(yself, j, xc);
00791         }
00792         xcprev = xc;
00793     }
00794 
00795     /* assert( i == 2 * hsize ); */
00796 
00797     if (xcnext >= xcprev) {
00798         /* Last data point is an extreme */
00799         /* FIXME: This could also be an error, recoverable by caller by
00800            increasing hsize */
00801         int j;
00802 
00803         imax++;
00804 
00805         if (cpl_bivector_get_size(self) < imax) {
00806             cpl_vector_set_size(xself, imax);
00807             cpl_vector_set_size(yself, imax);
00808         }
00809 
00810         for (j = imax-1; j > 0; j--) {
00811             if (xcnext <= cpl_vector_get(yself, j-1)) break;
00812             cpl_vector_set(xself, j, cpl_vector_get(xself, j-1));
00813             cpl_vector_set(yself, j, cpl_vector_get(yself, j-1));
00814         }
00815         cpl_vector_set(xself, j, i - hsize);
00816         cpl_vector_set(yself, j, xcnext);
00817 
00818     }
00819 
00820     if (doplot) {
00821         cpl_vector * xvals = cpl_vector_new(1 + 2 * hsize);
00822         cpl_bivector * bcorr = cpl_bivector_wrap_vectors(xvals, xcorr);
00823         double x = (double)-hsize;
00824         char * title = cpl_sprintf("t 'Cross-correlation of shifted %d-pixel "
00825                                    "spectrum (XCmax=%g at %d)' w linespoints",
00826                                    nobs, cpl_vector_get(xcorr, ixc),
00827                                    ixc - hsize);
00828 
00829         for (i = 0; i < 1 + 2 * hsize; i++, x += 1.0) {
00830             cpl_vector_set(xvals, i, x);
00831         }
00832 
00833         cpl_plot_bivector("set grid;set xlabel 'Offset [pixel]';", title,
00834                              "", bcorr);
00835         cpl_bivector_unwrap_vectors(bcorr);
00836         cpl_vector_delete(xvals);
00837         cpl_free(title);
00838     }
00839 
00840     if (pxc != NULL) *pxc = cpl_vector_get(xcorr, hsize);
00841 
00842     cpl_vector_delete(xcorr);
00843 
00844     if (imax < 1) {
00845         error = CPL_ERROR_DATA_NOT_FOUND;
00846     } else if (cpl_bivector_get_size(self) > imax) {
00847         cpl_vector_set_size(xself, imax);
00848         cpl_vector_set_size(yself, imax);
00849     }
00850 
00851     /* Propagate error, if any */
00852     return cpl_error_set(cpl_func, error);
00853 }
00854 
00855 /*----------------------------------------------------------------------------*/
00868 /*----------------------------------------------------------------------------*/
00869 cpl_error_code
00870 irplib_polynomial_shift_1d_from_correlation(cpl_polynomial * self,
00871                                             const cpl_vector * obs,
00872                                             irplib_base_spectrum_model * model,
00873                                             cpl_error_code (*filler)
00874                                             (cpl_vector *,
00875                                              const cpl_polynomial *,
00876                                              irplib_base_spectrum_model *, int),
00877                                             int hsize,
00878                                             cpl_boolean doplot,
00879                                             double * pxc)
00880 {
00881 
00882     const int      nobs   = cpl_vector_get_size(obs);
00883     const int      nmodel = 2 * hsize + nobs;
00884     cpl_vector   * mspec1d;
00885     cpl_vector   * xcorr;
00886     cpl_error_code error;
00887     int            ixc, xxc;
00888     double         xc;
00889 
00890     cpl_ensure_code(self   != NULL, CPL_ERROR_NULL_INPUT);
00891     cpl_ensure_code(obs    != NULL, CPL_ERROR_NULL_INPUT);
00892     cpl_ensure_code(model  != NULL, CPL_ERROR_NULL_INPUT);
00893     cpl_ensure_code(filler != NULL, CPL_ERROR_NULL_INPUT);
00894     cpl_ensure_code(hsize  >  0,    CPL_ERROR_ILLEGAL_INPUT);
00895 
00896     mspec1d = cpl_vector_new(nmodel);
00897 
00898     if (filler(mspec1d, self, model, hsize)) {
00899         cpl_vector_delete(mspec1d);
00900         cpl_ensure_code(0, cpl_error_get_code());
00901     }
00902 
00903     /* Should not be able to fail now */
00904     xcorr = cpl_vector_new(1 + 2 * hsize);
00905     ixc = cpl_vector_correlate(xcorr, mspec1d, obs);
00906 
00907 #ifdef IRPLIB_SPC_DUMP
00908     /* Need irplib_wavecal.c rev. 1.12 through 1.15 */
00909     irplib_polynomial_dump_corr_step(self, xcorr, "Shift");
00910 #endif
00911 
00912     cpl_vector_delete(mspec1d);
00913 
00914     xxc = ixc - hsize;
00915 
00916     error = cpl_polynomial_shift_1d(self, 0, (double)xxc);
00917 
00918     xc = cpl_vector_get(xcorr, ixc);
00919 
00920     cpl_msg_info(cpl_func, "Shifting %d pixels (%g < %g)", xxc,
00921                  cpl_vector_get(xcorr, hsize), xc);
00922 
00923     if (doplot) {
00924         cpl_vector * xvals = cpl_vector_new(1 + 2 * hsize);
00925         cpl_bivector * bcorr = cpl_bivector_wrap_vectors(xvals, xcorr);
00926         int i;
00927         double x = (double)-hsize;
00928         char * title = cpl_sprintf("t 'Cross-correlation of shifted %d-pixel "
00929                                    "spectrum (XCmax=%g at %d)' w linespoints",
00930                                    nobs, cpl_vector_get(xcorr, ixc), xxc);
00931 
00932         for (i = 0; i < 1 + 2 * hsize; i++, x += 1.0) {
00933             cpl_vector_set(xvals, i, x);
00934         }
00935 
00936         cpl_plot_bivector("set grid;set xlabel 'Offset [pixel]';", title,
00937                              "", bcorr);
00938         cpl_bivector_unwrap_vectors(bcorr);
00939         cpl_vector_delete(xvals);
00940         cpl_free(title);
00941     }
00942 
00943     cpl_vector_delete(xcorr);
00944 
00945     cpl_ensure_code(!error, error);
00946 
00947     if (pxc != NULL) *pxc = xc;
00948 
00949     return CPL_ERROR_NONE;
00950 
00951 }
00952 
00953 
00954 /*----------------------------------------------------------------------------*/
00974 /*----------------------------------------------------------------------------*/
00975 cpl_error_code
00976 irplib_vector_fill_line_spectrum_model(cpl_vector * self,
00977                                        cpl_vector * linepix,
00978                                        cpl_vector * erftmp,
00979                                        const cpl_polynomial * disp,
00980                                        const cpl_bivector * lines,
00981                                        double wslit,
00982                                        double wfwhm,
00983                                        double xtrunc,
00984                                        int hsize,
00985                                        cpl_boolean dofast,
00986                                        cpl_boolean dolog,
00987                                        unsigned * pulines)
00988 {
00989 
00990     cpl_errorstate     prestate;
00991     const double       sigma = wfwhm * CPL_MATH_SIG_FWHM;
00992     const cpl_vector * xlines  = cpl_bivector_get_x_const(lines);
00993     const double     * dxlines = cpl_vector_get_data_const(xlines);
00994     const double     * dylines = cpl_bivector_get_y_data_const(lines);
00995     double           * plinepix
00996         = linepix ? cpl_vector_get_data(linepix) : NULL;
00997     const int          nlines  = cpl_vector_get_size(xlines);
00998     const int          nself   = cpl_vector_get_size(self);
00999     double           * dself   = cpl_vector_get_data(self);
01000     cpl_polynomial   * dispi;
01001     double           * profile = NULL;
01002     const cpl_size     i0 = 0;
01003     const double       p0 = cpl_polynomial_get_coeff(disp, &i0);
01004     double             wl;
01005     double             xpos = (double)(1-hsize)-xtrunc;
01006     const double       xmax = (double)(nself-hsize)+xtrunc;
01007     double             xderiv, xextreme;
01008     cpl_error_code     error = CPL_ERROR_NONE;
01009     int                iline;
01010     unsigned           ulines = 0;
01011 
01012     cpl_ensure_code(self    != NULL, CPL_ERROR_NULL_INPUT);
01013     cpl_ensure_code(disp    != NULL, CPL_ERROR_NULL_INPUT);
01014     cpl_ensure_code(lines   != NULL, CPL_ERROR_NULL_INPUT);
01015 
01016     cpl_ensure_code(wslit  >  0.0, CPL_ERROR_ILLEGAL_INPUT);
01017     cpl_ensure_code(wfwhm  >  0.0, CPL_ERROR_ILLEGAL_INPUT);
01018     cpl_ensure_code(hsize  >= 0,   CPL_ERROR_ILLEGAL_INPUT);
01019     cpl_ensure_code(xtrunc >  0.0, CPL_ERROR_ILLEGAL_INPUT);
01020     cpl_ensure_code(nself  > 2 * hsize, CPL_ERROR_ILLEGAL_INPUT);
01021 
01022     cpl_ensure_code(cpl_polynomial_get_dimension(disp) == 1,
01023                     CPL_ERROR_ILLEGAL_INPUT);
01024     cpl_ensure_code(cpl_polynomial_get_degree(disp) > 0,
01025                     CPL_ERROR_ILLEGAL_INPUT);
01026 
01027     /* The smallest wavelength contributing to the spectrum. */
01028     wl = cpl_polynomial_eval_1d(disp, xpos, &xderiv);
01029 
01030     if (wl <= 0.0) return
01031         cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_INPUT, __FILE__,
01032                                     __LINE__, "Non-positive wavelength at x=%g: "
01033                                     "P(x)=%g, P'(x)=%g", xpos, wl, xderiv);
01034 
01035     if (xderiv <= 0.0) return
01036         cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_INPUT, __FILE__,
01037                                     __LINE__, "Non-increasing dispersion at "
01038                                     "x=%g: P'(x)=%g, P(x)=%g", xpos, xderiv, wl);
01039 
01040     /* Find the 1st line */
01041     iline = cpl_vector_find(xlines, wl);
01042 
01043     /* The first line must be at least at wl */
01044     if (dxlines[iline] < wl) iline++;
01045 
01046     if (iline >= nlines) return
01047         cpl_error_set_message_macro(cpl_func, CPL_ERROR_DATA_NOT_FOUND, __FILE__,
01048                                     __LINE__, "The %d-line catalogue has only "
01049                                     "lines below P(%g)=%g > %g", nlines, xpos,
01050                                     wl, dxlines[nlines-1]);
01051 
01052     memset(dself, 0, nself * sizeof(double));
01053 
01054     dispi = cpl_polynomial_duplicate(disp);
01055 
01056     /* Verify monotony of dispersion */
01057     cpl_polynomial_derivative(dispi, 0);
01058 
01059     prestate = cpl_errorstate_get();
01060 
01061     if (cpl_polynomial_solve_1d(dispi, 0.5*(nlines+1), &xextreme, 1)) {
01062         cpl_errorstate_set(prestate);
01063     } else if (xpos < xextreme && xextreme < xmax) {
01064         cpl_polynomial_delete(dispi);
01065         return cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
01066                                            __FILE__, __LINE__, "Non-monotone "
01067                                            "dispersion at x=%g: P'(x)=0, "
01068                                            "P(x)=%g", xextreme,
01069                                            cpl_polynomial_eval_1d(disp, xextreme,
01070                                                                   NULL));
01071     }
01072 
01073     if (dofast) {
01074         const int npix = 1+(int)xtrunc;
01075 
01076         if (erftmp != NULL && cpl_vector_get_size(erftmp) == npix &&
01077             cpl_vector_get(erftmp, 0) > 0.0) {
01078             profile = cpl_vector_get_data(erftmp);
01079         } else {
01080 
01081             const double yval =  0.5 / wslit;
01082             const double x0p  =  0.5 * wslit + 0.5;
01083             const double x0n  = -0.5 * wslit + 0.5;
01084             double       x1diff
01085                 = irplib_erf_antideriv(x0p, sigma)
01086                 - irplib_erf_antideriv(x0n, sigma);
01087             int          ipix;
01088 
01089             if (erftmp == NULL) {
01090                 profile = (double*)cpl_malloc(sizeof(double)*(size_t)npix);
01091             } else {
01092                 cpl_vector_set_size(erftmp, npix);
01093                 profile = cpl_vector_get_data(erftmp);
01094             }
01095 
01096             profile[0] = 2.0 * yval * x1diff;
01097 
01098             for (ipix = 1; ipix < npix; ipix++) {
01099                 const double x1 = (double)ipix;
01100                 const double x1p = x1 + 0.5 * wslit + 0.5;
01101                 const double x1n = x1 - 0.5 * wslit + 0.5;
01102                 const double x0diff = x1diff;
01103 
01104                 x1diff = irplib_erf_antideriv(x1p, sigma)
01105                     - irplib_erf_antideriv(x1n, sigma);
01106 
01107                 profile[ipix] = yval * (x1diff - x0diff);
01108 
01109             }
01110         }
01111     }
01112 
01113     cpl_polynomial_copy(dispi, disp);
01114 
01115     /* FIXME: A custom version of cpl_polynomial_solve_1d() which returns
01116        P'(xpos) can be used for the 1st NR-iteration. */
01117     /* Further, the sign of P'(xpos) could be checked for all lines. */
01118     /* Perform 1st NR-iteration in solving for P(xpos) = dxlines[iline] */
01119     xpos -= (wl - dxlines[iline]) / xderiv;
01120 
01121     /* Iterate through the lines */
01122     for (; !error && iline < nlines; iline++) {
01123 
01124         /* Lines may have a non-physical intensity (e.g. zero) to indicate some
01125            property of the line, e.g. unknown intensity due to blending */
01126         if (dylines[iline] <= 0.0) continue;
01127 
01128         /* Use 1st guess, if available (Use 0.0 to flag unavailable) */
01129         if (plinepix != NULL && plinepix[iline] > 0.0) xpos = plinepix[iline];
01130 
01131         if (xpos > xmax) xpos = xmax; /* FIXME: Better to limit xpos ? */
01132 
01133         /* Find the (sub-) pixel position of the line */
01134         error = cpl_polynomial_set_coeff(dispi, &i0, p0 - dxlines[iline]) ||
01135             cpl_polynomial_solve_1d(dispi, xpos, &xpos, 1);
01136 
01137         if (xpos > xmax) {
01138             if (error) {
01139                 error = 0;
01140                 cpl_msg_debug(cpl_func, "Stopping spectrum fill at line %d/%d "
01141                              "at xpos=%g > xmax=%g",
01142                              iline, nlines, xpos, xmax);
01143                 cpl_errorstate_dump(prestate, CPL_FALSE,
01144                                     irplib_errorstate_dump_debug);
01145                 cpl_errorstate_set(prestate);
01146             }
01147             break;
01148         } else if (error) {
01149             if (linepix != NULL && ulines) (void)cpl_vector_fill(linepix, 0.0);
01150             (void)cpl_error_set_message_macro(cpl_func, cpl_error_get_code(),
01151                                               __FILE__, __LINE__,
01152                                               "Could not find pixel-position "
01153                                               "of line %d/%d at wavelength=%g."
01154                                               " xpos=%g, xmax=%g",
01155                                               iline, nlines, dxlines[iline],
01156                                               xpos, xmax);
01157             break;
01158         } else if (dofast) {
01159             const double frac  = fabs(xpos - floor(xpos));
01160 #ifdef IRPLIB_WAVECAL_FAST_FAST
01161             const double frac0 = 1.0 - frac; /* Weight opposite of distance */
01162 #else
01163             /* Center intensity correctly */
01164             const double ep1pw = irplib_erf_antideriv(frac + 0.5 * wslit, sigma);
01165             const double en1pw = irplib_erf_antideriv(frac + 0.5 * wslit - 1.0,
01166                                                       sigma);
01167             const double ep1nw = irplib_erf_antideriv(frac - 0.5 * wslit, sigma);
01168             const double en1nw = irplib_erf_antideriv(frac - 0.5 * wslit - 1.0,
01169                                                       sigma);
01170             const double frac0
01171                 = (en1nw - en1pw) / (ep1pw - en1pw - ep1nw + en1nw);
01172 
01173 #endif
01174             const double frac1 = 1.0 - frac0;
01175             const double yval0 = frac0 * dylines[iline];
01176             const double yval1 = frac1 * dylines[iline];
01177             const int    npix  = 1+(int)xtrunc;
01178             int          ipix;
01179             int          i0n    = hsize - 1 + floor(xpos);
01180             int          i0p    = i0n;
01181             int          i1n    = i0n + 1;
01182             int          i1p    = i1n;
01183             cpl_boolean  didline = CPL_FALSE;
01184 
01185 
01186             /* Update 1st guess for next time, if available */
01187             if (plinepix != NULL) plinepix[iline] =  xpos;
01188 
01189             if (frac0 < 0.0) {
01190                 (void)cpl_error_set_message_macro(cpl_func,
01191                                                   CPL_ERROR_UNSPECIFIED,
01192                                                   __FILE__, __LINE__,
01193                                                   "Illegal split at x=%g: %g + "
01194                                                   "%g = 1", xpos, frac0, frac1);
01195 #ifdef IRPLIB_WAVEVAL_DEBUG
01196             } else {
01197                 cpl_msg_warning(cpl_func,"profile split at x=%g: %g + %g = 1",
01198                                 xpos, frac0, frac1);
01199 #endif
01200             }
01201 
01202             for (ipix = 0; ipix < npix; ipix++, i0n--, i0p++, i1n--, i1p++) {
01203 
01204                 if (i0n >= 0 && i0n < nself) {
01205                     dself[i0n] += yval0 * profile[ipix];
01206                     didline = CPL_TRUE;
01207                 }
01208                 if (i1n >= 0 && i1n < nself && ipix + 1 < npix) {
01209                     dself[i1n] += yval1 * profile[ipix+1];
01210                     didline = CPL_TRUE;
01211                 }
01212 
01213                 if (ipix == 0) continue;
01214 
01215                 if (i0p >= 0 && i0p < nself) {
01216                     dself[i0p] += yval0 * profile[ipix];
01217                     didline = CPL_TRUE;
01218                 }
01219                 if (i1p >= 0 && i1p < nself && ipix + 1 < npix) {
01220                     dself[i1p] += yval1 * profile[ipix+1];
01221                     didline = CPL_TRUE;
01222                 }
01223             }
01224 
01225             if (didline) ulines++;
01226 
01227         } else {
01228             const double yval = 0.5 * dylines[iline] / wslit;
01229             const int    ifirst = IRPLIB_MAX((int)(xpos-xtrunc+0.5), 1-hsize);
01230             const int    ilast  = IRPLIB_MIN((int)(xpos+xtrunc), nself-hsize);
01231             int          ipix;
01232             const double x0 = (double)ifirst - xpos;
01233             const double x0p = x0 + 0.5*wslit - 0.5;
01234             const double x0n = x0 - 0.5*wslit - 0.5;
01235             double       x1diff
01236                 = irplib_erf_antideriv(x0p, sigma)
01237                 - irplib_erf_antideriv(x0n, sigma);
01238 
01239             /* Update 1st guess for next time, if available */
01240             if (plinepix != NULL) plinepix[iline] =  xpos;
01241 
01242             if (ilast >= ifirst) ulines++;
01243 
01244             for (ipix = ifirst; ipix <= ilast; ipix++) {
01245                 const double x1 = (double)ipix - xpos;
01246                 const double x1p = x1 + 0.5*wslit + 0.5;
01247                 const double x1n = x1 - 0.5*wslit + 0.5;
01248                 const double x0diff = x1diff;
01249 
01250                 x1diff = irplib_erf_antideriv(x1p, sigma)
01251                     - irplib_erf_antideriv(x1n, sigma);
01252 
01253                 dself[ipix+hsize-1] += yval * (x1diff - x0diff);
01254 
01255             }
01256         }
01257     }
01258 
01259     cpl_polynomial_delete(dispi);
01260     if (erftmp == NULL) cpl_free(profile);
01261 
01262     cpl_ensure_code(!error, cpl_error_get_code());
01263 
01264     if (dolog) {
01265         int i;
01266         for (i = 0; i < nself; i++) {
01267             dself[i] = dself[i] > 0.0 ? log(1.0 + dself[i]) : 0.0;
01268         }
01269     }
01270 
01271     if (!ulines) return
01272         cpl_error_set_message_macro(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
01273                                     __FILE__, __LINE__, "The %d-line "
01274                                     "catalogue has no lines in the range "
01275                                     "%g -> P(%g)=%g", nlines, wl, xmax,
01276                                     cpl_polynomial_eval_1d(disp, xmax, NULL));
01277 
01278     if (pulines != NULL) *pulines = ulines;
01279 
01280     return CPL_ERROR_NONE;
01281 }
01282 
01283 /*----------------------------------------------------------------------------*/
01292 /*----------------------------------------------------------------------------*/
01293 inline double irplib_erf_antideriv(double x, double sigma)
01294 {
01295     return x * erf( x / (sigma * CPL_MATH_SQRT2))
01296        + 2.0 * sigma/CPL_MATH_SQRT2PI * exp(-0.5 * x * x / (sigma * sigma));
01297 }
01298 
01299 
01300 #ifdef HAVE_GSL
01301 
01302 /*----------------------------------------------------------------------------*/
01310 /*----------------------------------------------------------------------------*/
01311 static double irplib_gsl_correlation(const gsl_vector * self, void * data)
01312 {
01313 
01314     irplib_multimin * mindata = (irplib_multimin *)data;
01315     cpl_errorstate prestate = cpl_errorstate_get();
01316     int nobs, nmodel, ndiff;
01317     cpl_size i;
01318 
01319     cpl_ensure(self != NULL, CPL_ERROR_NULL_INPUT, GSL_NAN);
01320     cpl_ensure(data != NULL, CPL_ERROR_NULL_INPUT, GSL_NAN);
01321 
01322     cpl_ensure(mindata->filler   != NULL, CPL_ERROR_NULL_INPUT, GSL_NAN);
01323     cpl_ensure(mindata->observed != NULL, CPL_ERROR_NULL_INPUT, GSL_NAN);
01324     cpl_ensure(mindata->spectrum != NULL, CPL_ERROR_NULL_INPUT, GSL_NAN);
01325 
01326     nobs   = cpl_vector_get_size(mindata->observed);
01327     nmodel = cpl_vector_get_size(mindata->spectrum);
01328     ndiff  = nmodel - nobs;
01329 
01330     cpl_ensure((ndiff & 1) == 0, CPL_ERROR_ILLEGAL_INPUT, GSL_NAN);
01331 
01332     cpl_ensure(cpl_vector_get_size(mindata->vxc) == 1 + ndiff,
01333                CPL_ERROR_ILLEGAL_INPUT, GSL_NAN);
01334 
01335     ndiff /= 2;
01336 
01337     for (i=0; i < (cpl_size)self->size; i++) {
01338         const double value = gsl_vector_get(self, (size_t)i);
01339         cpl_polynomial_set_coeff(mindata->disp1d, &i, value);
01340     }
01341 
01342     if (mindata->filler(mindata->spectrum, mindata->disp1d,
01343                         mindata->param, ndiff)
01344         || !cpl_errorstate_is_equal(prestate)) {
01345 
01346         /* The fill failed. Ensure the discarding of this candidate by
01347            setting the cross-correlation to its minimum possible value. */
01348 
01349         (void)cpl_vector_fill(mindata->vxc, -1.0);
01350 
01351         mindata->maxxc = ndiff;
01352 
01353         if (!cpl_errorstate_is_equal(prestate)) {
01354                 cpl_msg_debug(cpl_func, "Spectrum fill failed:");
01355                 cpl_errorstate_dump(prestate, CPL_FALSE,
01356                                     irplib_errorstate_dump_debug);
01357                 cpl_errorstate_set(prestate);
01358         }
01359     } else {
01360 
01361         mindata->maxxc = cpl_vector_correlate(mindata->vxc,
01362                                               mindata->spectrum,
01363                                               mindata->observed);
01364     }
01365 
01366 #ifdef IRPLIB_SPC_DUMP
01367     /* Need irplib_wavecal.c rev. 1.12 through 1.15 */
01368     irplib_polynomial_dump_corr_step(mindata->disp1d, mindata->vxc,
01369                                      "Optimize");
01370 #endif
01371 
01372     mindata->xc = cpl_vector_get(mindata->vxc, ndiff);
01373 
01374     if (mindata->maxxc != ndiff &&
01375         cpl_vector_get(mindata->vxc, mindata->maxxc) > mindata->mxc) {
01376         const irplib_base_spectrum_model * arclamp
01377             = (const irplib_base_spectrum_model *)mindata->param;
01378 
01379         if (mindata->mdisp == NULL) {
01380             mindata->mdisp = cpl_polynomial_duplicate(mindata->disp1d);
01381         } else {
01382             cpl_polynomial_copy(mindata->mdisp, mindata->disp1d);
01383         }
01384         mindata->mxc = cpl_vector_get(mindata->vxc, mindata->maxxc);
01385         mindata->ishift = mindata->maxxc - ndiff;
01386         cpl_msg_debug(cpl_func, "Local maximum: %g(%d) > %g(%d) (cost=%u:%u. "
01387                       "lines=%u)", mindata->mxc, mindata->maxxc, mindata->xc,
01388                       ndiff, arclamp->cost, arclamp->xcost, arclamp->ulines);
01389     }
01390 
01391     return -mindata->xc;
01392 }
01393 
01394 #endif
01395 
01396 /*----------------------------------------------------------------------------*/
01419 /*----------------------------------------------------------------------------*/
01420 cpl_error_code
01421 irplib_polynomial_find_1d_from_correlation_all(cpl_polynomial * self,
01422                                                int maxdeg,
01423                                                const cpl_vector * obs,
01424                                                int nmaxima,
01425                                                int linelim,
01426                                                irplib_base_spectrum_model* model,
01427                                                cpl_error_code (* filler)
01428                                                (cpl_vector *,
01429                                                 const cpl_polynomial *,
01430                                                 irplib_base_spectrum_model *,
01431                                                 int),
01432                                                double pixtol,
01433                                                double pixstep,
01434                                                int hsize,
01435                                                int maxite,
01436                                                int maxfail,
01437                                                int maxcont,
01438                                                cpl_boolean doplot,
01439                                                double * pxc)
01440 {
01441 
01442 #ifdef HAVE_GSL
01443 
01444     cpl_errorstate     prestate = cpl_errorstate_get();
01445     cpl_polynomial   * start;
01446     cpl_polynomial   * cand;
01447     cpl_polynomial   * backup;
01448     cpl_error_code     error = CPL_ERROR_NONE;
01449     double             xc;
01450     cpl_bivector     * xtshift = cpl_bivector_new(nmaxima ? nmaxima : 1);
01451     const cpl_vector * xtshiftx = cpl_bivector_get_x_const(xtshift);
01452     const cpl_vector * xtshifty = cpl_bivector_get_y_const(xtshift);
01453     int                nshift;
01454     int                imaximum = -1;
01455     int                imaxima;
01456 
01457 #endif
01458 
01459     cpl_ensure_code(self   != NULL, CPL_ERROR_NULL_INPUT);
01460     cpl_ensure_code(obs    != NULL, CPL_ERROR_NULL_INPUT);
01461     cpl_ensure_code(model  != NULL, CPL_ERROR_NULL_INPUT);
01462     cpl_ensure_code(filler != NULL, CPL_ERROR_NULL_INPUT);
01463     cpl_ensure_code(pxc    != NULL, CPL_ERROR_NULL_INPUT);
01464 
01465     cpl_ensure_code(cpl_polynomial_get_dimension(self) == 1,
01466                     CPL_ERROR_ILLEGAL_INPUT);
01467 
01468     cpl_ensure_code(cpl_polynomial_get_degree(self) > 0,
01469                     CPL_ERROR_ILLEGAL_INPUT);
01470 
01471     cpl_ensure_code(maxdeg  >=  0, CPL_ERROR_ILLEGAL_INPUT);
01472     cpl_ensure_code(pixtol  > 0.0, CPL_ERROR_ILLEGAL_INPUT);
01473     cpl_ensure_code(pixstep > 0.0, CPL_ERROR_ILLEGAL_INPUT);
01474     cpl_ensure_code(hsize   >=  0, CPL_ERROR_ILLEGAL_INPUT);
01475     cpl_ensure_code(maxite  >=  0, CPL_ERROR_ILLEGAL_INPUT);
01476     cpl_ensure_code(nmaxima >=  0, CPL_ERROR_ILLEGAL_INPUT);
01477     cpl_ensure_code(maxfail >   0, CPL_ERROR_ILLEGAL_INPUT);
01478     cpl_ensure_code(maxcont >   0, CPL_ERROR_ILLEGAL_INPUT);
01479     cpl_ensure_code(linelim >=  0, CPL_ERROR_ILLEGAL_INPUT);
01480 
01481 #ifndef HAVE_GSL
01482     /* Avoid unused variable warning */
01483     cpl_ensure_code(doplot == CPL_TRUE || doplot == CPL_FALSE,
01484                     CPL_ERROR_ILLEGAL_INPUT);
01485     return cpl_error_set_message(cpl_func, CPL_ERROR_UNSUPPORTED_MODE,
01486                                  "GSL is not available");
01487 #else
01488 
01489     if (irplib_bivector_find_shift_from_correlation(xtshift, self, obs,
01490                                                     model, filler,
01491                                                     hsize, doplot, &xc)) {
01492         cpl_bivector_delete(xtshift);
01493         return cpl_error_set_where(cpl_func);
01494     }
01495 
01496     if (model->ulines > (unsigned)linelim) {
01497         /* The initial, optimal (integer) shift */
01498         const double xxc = cpl_vector_get(xtshiftx, 0);
01499         const double xc0 = cpl_vector_get(xtshifty, 0);
01500 
01501         cpl_msg_warning(cpl_func, "Doing only shift=%g pixels with lines=%u > "
01502                         "%d and XC=%g", xxc, model->ulines, linelim, xc0);
01503 
01504         cpl_polynomial_shift_1d(self, 0, xxc);
01505 
01506         *pxc = xc0;
01507 
01508         cpl_bivector_delete(xtshift);
01509 
01510         return CPL_ERROR_NONE;
01511     }
01512 
01513     start  = cpl_polynomial_duplicate(self);
01514     cand   = cpl_polynomial_new(1);
01515     backup = cpl_polynomial_new(1);
01516 
01517     /* Number of (local) maxima to use as starting point of the optimization */
01518     nshift = cpl_bivector_get_size(xtshift);
01519     if (nmaxima == 0 || nmaxima > nshift) nmaxima = nshift;
01520 
01521     cpl_msg_info(cpl_func, "Optimizing %d/%d local shift-maxima "
01522                  "(no-shift xc=%g. linelim=%d)", nmaxima, nshift, xc, linelim);
01523     if (cpl_msg_get_level() <= CPL_MSG_DEBUG)
01524         cpl_bivector_dump(xtshift, stdout);
01525 
01526     for (imaxima = 0; imaxima < nmaxima; imaxima++) {
01527         /* The initial, optimal (integer) shift */
01528         const double xxc = cpl_vector_get(xtshiftx, imaxima);
01529         double xtpixstep = pixstep;
01530         double xtpixtol  = pixtol;
01531         double xtxc;
01532         cpl_boolean ok = CPL_FALSE;
01533         int nfail;
01534 
01535 
01536         cpl_polynomial_copy(cand, start);
01537         cpl_polynomial_shift_1d(cand, 0, xxc);
01538         cpl_polynomial_copy(backup, cand);
01539 
01540         /* Increase tolerance until convergence */
01541         for (nfail = 0; nfail < maxfail; nfail++, xtpixtol *= 2.0,
01542                  xtpixstep *= 2.0) {
01543             int restart = maxcont;
01544             do {
01545                 error = irplib_polynomial_find_1d_from_correlation
01546                     (cand, maxdeg, obs, model,
01547                      filler, xtpixtol, xtpixstep, 2,
01548                      maxite, &xtxc);
01549             } while (error == CPL_ERROR_CONTINUE && --restart);
01550 
01551             if (!error) {
01552                 cpl_msg_debug(cpl_func, "XC(imax=%d/%d:xtpixtol=%g): %g "
01553                               "(cost=%u:%u)", 1+imaxima, nmaxima, xtpixtol,
01554                               xtxc, model->cost, model->xcost);
01555                 break;
01556             }
01557             cpl_msg_warning(cpl_func, "Increasing xtpixtol from %g (%g, imax="
01558                             "%d/%d)", xtpixtol, xtpixstep, 1+imaxima, nmaxima);
01559             if (model->ulines > (unsigned)linelim) {
01560                 cpl_msg_warning(cpl_func, "Stopping search-refinement via "
01561                                 "catalogue with %u lines > %u", model->ulines,
01562                                 linelim);
01563                 break;
01564             }
01565             cpl_polynomial_copy(cand, start);
01566         }
01567 
01568         /* Decrease tolerance until inconvergence, keep previous */
01569         for (; !error && xtpixtol > 0.0; xtpixtol *= 0.25, xtpixstep *= 0.5) {
01570             int restart = maxcont;
01571 
01572             cpl_polynomial_copy(backup, cand);
01573             do {
01574                 error = irplib_polynomial_find_1d_from_correlation
01575                     (cand, maxdeg, obs, model, filler,
01576                      xtpixtol, xtpixstep, 2, maxite, &xtxc);
01577             } while (error == CPL_ERROR_CONTINUE && --restart);
01578             if (error) break;
01579             ok = CPL_TRUE;
01580             cpl_msg_debug(cpl_func, "XC(imax=%d/%d:xtpixtol=%g): %g (cost=%u:%u"
01581                           ". ulines=%u)", 1+imaxima, nmaxima, xtpixtol, xtxc,
01582                           model->cost, model->xcost, model->ulines);
01583             if (model->ulines > (unsigned)linelim) {
01584                 cpl_msg_info(cpl_func, "Stopping search-refinement via "
01585                              "catalogue with %u lines > %u", model->ulines,
01586                              linelim);
01587                 break;
01588             }
01589         }
01590 
01591         if (error) {
01592             error = 0;
01593             cpl_errorstate_dump(prestate, CPL_FALSE,
01594                                 irplib_errorstate_dump_debug);
01595             cpl_errorstate_set(prestate);
01596             cpl_polynomial_copy(cand, backup);
01597         }
01598         if (ok && xtxc > xc) {
01599             imaximum = imaxima;
01600             cpl_polynomial_copy(self, cand);
01601             xc = xtxc;
01602 
01603             cpl_msg_info(cpl_func, "XC(imax=%d/%d): %g -> %g (initial-shift=%g. "
01604                          "cost=%u:%u. lines=%u)", 1+imaxima, nmaxima,
01605                          cpl_vector_get(xtshifty, imaxima), xtxc,
01606                          cpl_vector_get(xtshiftx, imaxima), model->cost,
01607                          model->xcost, model->ulines);
01608         } else {
01609             cpl_msg_info(cpl_func, "xc(imax=%d/%d): %g -> %g (initial-shift=%g. "
01610                          "cost=%u:%u. lines=%u)", 1+imaxima, nmaxima,
01611                          cpl_vector_get(xtshifty, imaxima), xtxc,
01612                          cpl_vector_get(xtshiftx, imaxima),
01613                          model->cost, model->xcost, model->ulines);
01614         }
01615     }
01616 
01617     cpl_polynomial_delete(start);
01618     cpl_polynomial_delete(backup);
01619     cpl_polynomial_delete(cand);
01620 
01621     if (imaximum < 0) {
01622         error = cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
01623                                       "XC could not be optimized over %d "
01624                                       "local shift-maxima (xc=%g)", nmaxima, xc);
01625     } else {
01626         cpl_msg_info(cpl_func, "Maximal XC=%g (up from %g, with initial pixel-"
01627                      "shift of %g) at %d/%d local shift-maximi", xc,
01628                      cpl_vector_get(xtshifty, imaximum),
01629                      cpl_vector_get(xtshiftx, imaximum),
01630                      1+imaximum, nmaxima);
01631 
01632         if (doplot) {
01633             irplib_plot_spectrum_and_model(obs, self, model, filler);
01634         }
01635 
01636         *pxc = xc;
01637     }
01638 
01639     cpl_bivector_delete(xtshift);
01640 
01641     return error;
01642 
01643 #endif
01644 
01645 }

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