28 #include "irplib_wavecal_impl.h"
31 #include "irplib_utils.h"
39 #include <gsl/gsl_multimin.h>
48 const cpl_vector * observed;
49 cpl_polynomial * disp1d;
50 cpl_vector * spectrum;
51 irplib_base_spectrum_model * param;
52 cpl_error_code (* filler)(cpl_vector *,
const cpl_polynomial *,
53 irplib_base_spectrum_model *);
58 cpl_polynomial * mdisp;
71 #define IRPLIB_MAX(A,B) ((A) > (B) ? (A) : (B))
72 #define IRPLIB_MIN(A,B) ((A) < (B) ? (A) : (B))
79 static double irplib_gsl_correlation(
const gsl_vector *,
void *);
83 irplib_polynomial_find_1d_from_correlation_(cpl_polynomial *,
int,
85 irplib_base_spectrum_model *,
88 const cpl_polynomial *,
89 irplib_base_spectrum_model *),
90 double,
double,
int,
int,
91 double *, cpl_boolean *);
117 const int nself = cpl_bivector_get_size(
self);
118 const double * px = cpl_bivector_get_x_data_const(
self);
119 const double * py = cpl_bivector_get_y_data_const(
self);
123 cpl_ensure(
self != NULL, CPL_ERROR_NULL_INPUT, -1);
124 cpl_ensure(x_min <= x_max, CPL_ERROR_ILLEGAL_INPUT, -2);
127 while (i < nself && px[i] < x_min) i++;
128 while (i < nself && px[i] < x_max)
129 if (py[i++] > 0) npos++;
146 const cpl_image * imgwave,
147 int fitdeg,
double * presid)
150 const int nx = cpl_image_get_size_x(imgwave);
151 const int ny = cpl_image_get_size_y(imgwave);
152 const int nbad = cpl_image_count_rejected(imgwave);
153 const int nsamp = nx * ny - nbad;
159 const cpl_size nfitdeg = (cpl_size)fitdeg;
163 cpl_ensure_code(
self != NULL, CPL_ERROR_NULL_INPUT);
164 cpl_ensure_code(imgwave != NULL, CPL_ERROR_NULL_INPUT);
165 cpl_ensure_code(presid != NULL, CPL_ERROR_NULL_INPUT);
166 cpl_ensure_code(fitdeg > 0, CPL_ERROR_ILLEGAL_INPUT);
168 cpl_ensure_code(cpl_polynomial_get_dimension(
self) == 2,
169 CPL_ERROR_ILLEGAL_INPUT);
171 xy_pos = cpl_matrix_new(2, nsamp);
172 xdata = cpl_matrix_get_data(xy_pos);
173 ydata = xdata + nsamp;
175 dwlen = (
double*)cpl_malloc(nsamp *
sizeof(
double));
176 wlen = cpl_vector_wrap(nsamp, dwlen);
178 for (i=1; i <= nx; i++) {
179 for (j=1; j <= ny; j++) {
181 const double value = cpl_image_get(imgwave, i, j, &is_bad);
191 cpl_msg_info(cpl_func,
"Fitting 2D polynomial to %d X %d image, ignoring "
192 "%d poorly calibrated pixels", nx, ny, nbad);
194 if (cpl_polynomial_fit(
self, xy_pos, NULL, wlen, NULL, CPL_FALSE, NULL,
195 &nfitdeg) == CPL_ERROR_NONE && presid != NULL) {
196 cpl_vector_fill_polynomial_fit_residual(wlen, wlen, NULL,
self, xy_pos,
198 *presid = cpl_vector_product(wlen, wlen)/nsamp;
200 cpl_matrix_delete(xy_pos);
201 cpl_vector_delete(wlen);
203 cpl_ensure_code(k == nsamp, CPL_ERROR_UNSPECIFIED);
205 return CPL_ERROR_NONE;
231 const cpl_vector * obs,
232 irplib_base_spectrum_model * model,
233 cpl_error_code (* filler)
235 const cpl_polynomial *,
236 irplib_base_spectrum_model *),
243 cpl_boolean restart = CPL_FALSE;
244 const cpl_error_code error = irplib_polynomial_find_1d_from_correlation_
245 (
self, maxdeg, obs, model, filler, pixtol, pixstep, hsize, maxite, pxc,
248 return error ? cpl_error_set_where(cpl_func) :
249 (restart ? cpl_error_set(cpl_func, CPL_ERROR_CONTINUE)
275 static cpl_error_code
276 irplib_polynomial_find_1d_from_correlation_(cpl_polynomial *
self,
278 const cpl_vector * obs,
279 irplib_base_spectrum_model * model,
280 cpl_error_code (* filler)
282 const cpl_polynomial *,
283 irplib_base_spectrum_model *),
289 cpl_boolean * prestart)
293 const gsl_multimin_fminimizer_type * T = gsl_multimin_fminimizer_nmsimplex;
294 gsl_multimin_fminimizer * minimizer;
295 gsl_multimin_function my_func;
296 irplib_multimin data;
297 gsl_vector * dispgsl;
298 gsl_vector * stepsize;
299 gsl_vector * dispprev;
300 int status = GSL_CONTINUE;
301 const int nobs = cpl_vector_get_size(obs);
302 const cpl_size nfit = maxdeg + 1;
303 cpl_errorstate prestate = cpl_errorstate_get();
305 const double wlstep =
306 cpl_polynomial_eval_1d_diff(
self, 0.5 * (nobs + pixstep),
307 0.5 * (nobs - pixstep), NULL);
308 double wlstepi = wlstep;
315 cpl_ensure_code(prestart != NULL, CPL_ERROR_NULL_INPUT);
316 *prestart = CPL_FALSE;
317 cpl_ensure_code(
self != NULL, CPL_ERROR_NULL_INPUT);
318 cpl_ensure_code(obs != NULL, CPL_ERROR_NULL_INPUT);
319 cpl_ensure_code(model != NULL, CPL_ERROR_NULL_INPUT);
320 cpl_ensure_code(filler != NULL, CPL_ERROR_NULL_INPUT);
321 cpl_ensure_code(pxc != NULL, CPL_ERROR_NULL_INPUT);
323 cpl_ensure_code(cpl_polynomial_get_dimension(
self) == 1,
324 CPL_ERROR_ILLEGAL_INPUT);
326 cpl_ensure_code(cpl_polynomial_get_degree(
self) > 0,
327 CPL_ERROR_ILLEGAL_INPUT);
329 cpl_ensure_code(maxdeg >= 0, CPL_ERROR_ILLEGAL_INPUT);
330 cpl_ensure_code(pixtol > 0.0, CPL_ERROR_ILLEGAL_INPUT);
331 cpl_ensure_code(pixstep > 0.0, CPL_ERROR_ILLEGAL_INPUT);
332 cpl_ensure_code(hsize >= 0, CPL_ERROR_ILLEGAL_INPUT);
333 cpl_ensure_code(maxite >= 0, CPL_ERROR_ILLEGAL_INPUT);
336 return cpl_error_set_message(cpl_func, CPL_ERROR_UNSUPPORTED_MODE,
337 "GSL is not available");
340 minimizer = gsl_multimin_fminimizer_alloc(T, (
size_t)nfit);
342 cpl_ensure_code(minimizer != NULL, CPL_ERROR_ILLEGAL_OUTPUT);
344 dispgsl = gsl_vector_alloc((
size_t)nfit);
345 stepsize = gsl_vector_alloc((
size_t)nfit);
346 dispprev = gsl_vector_alloc((
size_t)nfit);
348 for (i=0; i < nfit; i++) {
349 const double value = cpl_polynomial_get_coeff(
self, &i);
350 gsl_vector_set(dispgsl, (
size_t)i, value);
351 gsl_vector_set(stepsize, (
size_t)i, wlstepi);
352 wlstepi /= (double)nobs;
356 my_func.f = &irplib_gsl_correlation;
357 my_func.params = (
void *)(&data);
361 data.spectrum = cpl_vector_new(nobs + 2 * hsize);
362 data.vxc = cpl_vector_new(1 + 2 * hsize);
364 data.filler = filler;
370 gsl_multimin_fminimizer_set (minimizer, &my_func, dispgsl, stepsize);
372 for (iter = 0; status == GSL_CONTINUE && iter < maxite; iter++) {
374 const double fprev = minimizer->fval;
376 gsl_vector_memcpy(dispprev, minimizer->x);
377 status = gsl_multimin_fminimizer_iterate(minimizer);
379 if (status || !cpl_errorstate_is_equal(prestate))
break;
381 size = gsl_multimin_fminimizer_size (minimizer);
382 status = gsl_multimin_test_size (size, pixtol);
384 if (status == GSL_SUCCESS) {
385 cpl_msg_debug(cpl_func,
"converged to minimum at");
388 cpl_msg_debug(cpl_func,
"%5d %g df() = %g size = %g",
390 gsl_vector_get (minimizer->x, 0)
391 - gsl_vector_get (dispprev, 0),
392 minimizer->fval - fprev, size);
393 }
else if (nfit == 1) {
394 cpl_msg_debug(cpl_func,
"%5d %g %g df() = %g size = %g",
396 gsl_vector_get (minimizer->x, 0)
397 - gsl_vector_get (dispprev, 0),
398 gsl_vector_get (minimizer->x, 1)
399 - gsl_vector_get (dispprev, 1),
400 minimizer->fval - fprev, size);
402 cpl_msg_debug(cpl_func,
"%5d %g %g %g df() = %g size = %g",
404 gsl_vector_get (minimizer->x, 0)
405 - gsl_vector_get (dispprev, 0),
406 gsl_vector_get (minimizer->x, 1)
407 - gsl_vector_get (dispprev, 1),
408 gsl_vector_get (minimizer->x, 2)
409 - gsl_vector_get (dispprev, 2),
410 minimizer->fval - fprev, size);
415 if (status == GSL_SUCCESS && cpl_errorstate_is_equal(prestate)) {
416 if (data.mxc > -minimizer->fval) {
418 cpl_msg_warning(cpl_func,
"Local maximum: %g(%d) > %g",
419 data.mxc, data.ishift, -minimizer->fval);
420 cpl_polynomial_shift_1d(data.mdisp, 0, (
double)data.ishift);
421 cpl_polynomial_copy(
self, data.mdisp);
422 *prestart = CPL_TRUE;
424 *pxc = -minimizer->fval;
425 for (i=0; i < nfit; i++) {
426 const double value = gsl_vector_get(minimizer->x, i);
427 cpl_polynomial_set_coeff(
self, &i, value);
432 cpl_vector_delete(data.spectrum);
433 cpl_vector_delete(data.vxc);
434 cpl_polynomial_delete(data.mdisp);
435 gsl_multimin_fminimizer_free(minimizer);
436 gsl_vector_free(dispgsl);
437 gsl_vector_free(dispprev);
438 gsl_vector_free(stepsize);
440 cpl_ensure_code(status != GSL_CONTINUE, CPL_ERROR_CONTINUE);
441 cpl_ensure_code(status == GSL_SUCCESS, CPL_ERROR_DATA_NOT_FOUND);
442 cpl_ensure_code(cpl_errorstate_is_equal(prestate), cpl_error_get_code());
444 return CPL_ERROR_NONE;
480 const cpl_polynomial * disp,
481 irplib_base_spectrum_model * lsslamp)
484 irplib_line_spectrum_model * arclamp
485 = (irplib_line_spectrum_model *)lsslamp;
486 cpl_error_code error;
488 cpl_ensure_code(arclamp != NULL, CPL_ERROR_NULL_INPUT);
500 0, CPL_FALSE, CPL_FALSE,
502 cpl_ensure_code(!error, error);
506 return CPL_ERROR_NONE;
525 const cpl_polynomial * disp,
526 irplib_base_spectrum_model * lsslamp)
529 irplib_line_spectrum_model * arclamp
530 = (irplib_line_spectrum_model *)lsslamp;
531 cpl_error_code error;
533 cpl_ensure_code(arclamp != NULL, CPL_ERROR_NULL_INPUT);
545 0, CPL_FALSE, CPL_TRUE,
547 cpl_ensure_code(!error, error);
551 return CPL_ERROR_NONE;
571 const cpl_polynomial * disp,
572 irplib_base_spectrum_model * lsslamp)
575 irplib_line_spectrum_model * arclamp
576 = (irplib_line_spectrum_model *)lsslamp;
577 cpl_error_code error;
579 cpl_ensure_code(arclamp != NULL, CPL_ERROR_NULL_INPUT);
591 0, CPL_TRUE, CPL_FALSE,
593 cpl_ensure_code(!error, error);
597 return CPL_ERROR_NONE;
616 const cpl_polynomial * disp,
617 irplib_base_spectrum_model * lsslamp)
620 irplib_line_spectrum_model * arclamp
621 = (irplib_line_spectrum_model *)lsslamp;
622 cpl_error_code error;
624 cpl_ensure_code(arclamp != NULL, CPL_ERROR_NULL_INPUT);
636 0, CPL_TRUE, CPL_TRUE,
638 cpl_ensure_code(!error, error);
642 return CPL_ERROR_NONE;
658 const cpl_polynomial * disp1d,
659 irplib_base_spectrum_model * model,
660 cpl_error_code (* filler)
662 const cpl_polynomial *,
663 irplib_base_spectrum_model *))
666 cpl_errorstate prestate = cpl_errorstate_get();
668 cpl_vector * spectrum;
670 const int len = cpl_vector_get_size(
self);
675 cpl_ensure_code(
self != NULL, CPL_ERROR_NULL_INPUT);
676 cpl_ensure_code(disp1d != NULL, CPL_ERROR_NULL_INPUT);
677 cpl_ensure_code(model != NULL, CPL_ERROR_NULL_INPUT);
678 cpl_ensure_code(filler != NULL, CPL_ERROR_NULL_INPUT);
680 cpl_ensure_code(cpl_polynomial_get_dimension(disp1d) == 1,
681 CPL_ERROR_ILLEGAL_INPUT);
683 cpl_ensure_code(cpl_polynomial_get_degree(disp1d) > 0,
684 CPL_ERROR_ILLEGAL_INPUT);
686 wl = cpl_vector_new(len);
687 spectrum = cpl_vector_new(len);
688 vxc = cpl_vector_new(1);
690 error |= (int)cpl_vector_fill_polynomial(wl, disp1d, 1.0, 1.0);
691 error |= filler(spectrum, disp1d, model);
693 ixc = cpl_vector_correlate(vxc,
self, spectrum);
694 xc = cpl_vector_get(vxc, ixc);
696 maxval = cpl_vector_get_max(spectrum);
698 error |= cpl_vector_multiply_scalar(spectrum,
699 cpl_vector_get_max(
self)/maxval);
701 const cpl_vector * spair[] = {wl,
self, spectrum};
702 char * pre = cpl_sprintf(
"set grid;set xlabel 'Wavelength (%g -> %g)'; "
703 "set ylabel 'Intensity';", cpl_vector_get(wl, 0),
704 cpl_vector_get(wl, len-1));
705 char * title = cpl_sprintf(
"t 'Observed and modelled spectra (%d pixel "
706 "XC=%g) ' w linespoints", len, xc);
708 (void)cpl_plot_vectors(pre, title,
"", spair, 3);
713 cpl_vector_delete(wl);
714 cpl_vector_delete(spectrum);
715 cpl_vector_delete(vxc);
717 cpl_errorstate_set(prestate);
719 return CPL_ERROR_NONE;
745 const cpl_polynomial * disp,
746 const cpl_vector * obs,
747 irplib_base_spectrum_model * model,
748 cpl_error_code (*filler)
750 const cpl_polynomial *,
751 irplib_base_spectrum_model *),
757 const int nobs = cpl_vector_get_size(obs);
758 const int nmodel = 2 * hsize + nobs;
759 cpl_polynomial * shdisp;
760 cpl_vector * xself = cpl_bivector_get_x(
self);
761 cpl_vector * yself = cpl_bivector_get_y(
self);
762 cpl_vector * mspec1d;
764 cpl_error_code error = CPL_ERROR_NONE;
765 double xcprev, xcnext;
769 cpl_ensure_code(
self != NULL, CPL_ERROR_NULL_INPUT);
770 cpl_ensure_code(disp != NULL, CPL_ERROR_NULL_INPUT);
771 cpl_ensure_code(obs != NULL, CPL_ERROR_NULL_INPUT);
772 cpl_ensure_code(model != NULL, CPL_ERROR_NULL_INPUT);
773 cpl_ensure_code(filler != NULL, CPL_ERROR_NULL_INPUT);
774 cpl_ensure_code(hsize > 0, CPL_ERROR_ILLEGAL_INPUT);
776 shdisp = cpl_polynomial_duplicate(disp);
779 if (cpl_polynomial_shift_1d(shdisp, 0, -hsize)) {
780 cpl_polynomial_delete(shdisp);
781 return cpl_error_set_where(cpl_func);
784 mspec1d = cpl_vector_new(nmodel);
786 if (filler(mspec1d, shdisp, model)) {
787 cpl_vector_delete(mspec1d);
788 return cpl_error_set_where(cpl_func);
792 xcorr = cpl_vector_new(1 + 2 * hsize);
793 ixc = cpl_vector_correlate(xcorr, mspec1d, obs);
795 #ifdef IRPLIB_SPC_DUMP
797 irplib_polynomial_dump_corr_step(shdisp, xcorr,
"Shift");
800 cpl_vector_delete(mspec1d);
801 cpl_polynomial_delete(shdisp);
806 xcprev = cpl_vector_get(xcorr, i);
807 xcnext = cpl_vector_get(xcorr, i+1);
809 if (xcprev >= xcnext) {
815 cpl_vector_set(xself, 0, i - hsize);
816 cpl_vector_set(yself, 0, xcprev);
820 for (i = 1; i < 2 * hsize; i++) {
821 const double xc = xcnext;
822 xcnext = cpl_vector_get(xcorr, i+1);
823 if (xc >= xcprev && xc >= xcnext) {
829 if (cpl_bivector_get_size(
self) < imax) {
830 cpl_vector_set_size(xself, imax);
831 cpl_vector_set_size(yself, imax);
834 for (j = imax-1; j > 0; j--) {
835 if (xc <= cpl_vector_get(yself, j-1))
break;
836 cpl_vector_set(xself, j, cpl_vector_get(xself, j-1));
837 cpl_vector_set(yself, j, cpl_vector_get(yself, j-1));
839 cpl_vector_set(xself, j, i - hsize);
840 cpl_vector_set(yself, j, xc);
847 if (xcnext >= xcprev) {
855 if (cpl_bivector_get_size(
self) < imax) {
856 cpl_vector_set_size(xself, imax);
857 cpl_vector_set_size(yself, imax);
860 for (j = imax-1; j > 0; j--) {
861 if (xcnext <= cpl_vector_get(yself, j-1))
break;
862 cpl_vector_set(xself, j, cpl_vector_get(xself, j-1));
863 cpl_vector_set(yself, j, cpl_vector_get(yself, j-1));
865 cpl_vector_set(xself, j, i - hsize);
866 cpl_vector_set(yself, j, xcnext);
872 cpl_vector * xvals = cpl_vector_new(1 + 2 * hsize);
873 cpl_bivector * bcorr = cpl_bivector_wrap_vectors(xvals, xcorr);
874 double x = (double)-hsize;
875 char * title = cpl_sprintf(
"t 'Cross-correlation of shifted %d-pixel "
876 "spectrum (XCmax=%g at %d)' w linespoints",
877 nobs, cpl_vector_get(xcorr, ixc),
880 for (i = 0; i < 1 + 2 * hsize; i++, x += 1.0) {
881 cpl_vector_set(xvals, i, x);
884 cpl_plot_bivector(
"set grid;set xlabel 'Offset [pixel]';", title,
886 cpl_bivector_unwrap_vectors(bcorr);
887 cpl_vector_delete(xvals);
891 if (pxc != NULL) *pxc = cpl_vector_get(xcorr, hsize);
893 cpl_vector_delete(xcorr);
896 error = CPL_ERROR_DATA_NOT_FOUND;
897 }
else if (cpl_bivector_get_size(
self) > imax) {
898 cpl_vector_set_size(xself, imax);
899 cpl_vector_set_size(yself, imax);
903 return cpl_error_set(cpl_func, error);
922 const cpl_vector * obs,
923 irplib_base_spectrum_model * model,
924 cpl_error_code (*filler)
926 const cpl_polynomial *,
927 irplib_base_spectrum_model *),
933 const int nobs = cpl_vector_get_size(obs);
934 const int nmodel = 2 * hsize + nobs;
935 cpl_vector * mspec1d;
937 cpl_error_code error;
941 cpl_ensure_code(
self != NULL, CPL_ERROR_NULL_INPUT);
942 cpl_ensure_code(obs != NULL, CPL_ERROR_NULL_INPUT);
943 cpl_ensure_code(model != NULL, CPL_ERROR_NULL_INPUT);
944 cpl_ensure_code(filler != NULL, CPL_ERROR_NULL_INPUT);
945 cpl_ensure_code(hsize > 0, CPL_ERROR_ILLEGAL_INPUT);
948 cpl_ensure_code(!cpl_polynomial_shift_1d(
self, 0, -hsize),
949 cpl_error_get_code());
951 mspec1d = cpl_vector_new(nmodel);
953 if (filler(mspec1d,
self, model)) {
954 cpl_vector_delete(mspec1d);
955 cpl_ensure_code(0, cpl_error_get_code());
959 xcorr = cpl_vector_new(1 + 2 * hsize);
960 ixc = cpl_vector_correlate(xcorr, mspec1d, obs);
962 #ifdef IRPLIB_SPC_DUMP
964 irplib_polynomial_dump_corr_step(
self, xcorr,
"Shift");
967 cpl_vector_delete(mspec1d);
969 error = cpl_polynomial_shift_1d(
self, 0, (
double)ixc);
971 xc = cpl_vector_get(xcorr, ixc);
975 cpl_msg_info(cpl_func,
"Shifting %d pixels (%g < %g)", xxc,
976 cpl_vector_get(xcorr, hsize), xc);
979 cpl_vector * xvals = cpl_vector_new(1 + 2 * hsize);
980 cpl_bivector * bcorr = cpl_bivector_wrap_vectors(xvals, xcorr);
982 double x = (double)-hsize;
983 char * title = cpl_sprintf(
"t 'Cross-correlation of shifted %d-pixel "
984 "spectrum (XCmax=%g at %d)' w linespoints",
985 nobs, cpl_vector_get(xcorr, ixc), xxc);
987 for (i = 0; i < 1 + 2 * hsize; i++, x += 1.0) {
988 cpl_vector_set(xvals, i, x);
991 cpl_plot_bivector(
"set grid;set xlabel 'Offset [pixel]';", title,
993 cpl_bivector_unwrap_vectors(bcorr);
994 cpl_vector_delete(xvals);
998 cpl_vector_delete(xcorr);
1000 cpl_ensure_code(!error, error);
1002 if (pxc != NULL) *pxc = xc;
1004 return CPL_ERROR_NONE;
1032 cpl_vector * linepix,
1033 cpl_vector * erftmp,
1034 const cpl_polynomial * disp,
1035 const cpl_bivector * lines,
1045 cpl_errorstate prestate;
1046 const double sigma = wfwhm * CPL_MATH_SIG_FWHM;
1047 const cpl_vector * xlines = cpl_bivector_get_x_const(lines);
1048 const double * dxlines = cpl_vector_get_data_const(xlines);
1049 const double * dylines = cpl_bivector_get_y_data_const(lines);
1051 = linepix ? cpl_vector_get_data(linepix) : NULL;
1052 const int nlines = cpl_vector_get_size(xlines);
1053 const int nself = cpl_vector_get_size(
self);
1054 double * dself = cpl_vector_get_data(
self);
1055 cpl_polynomial * dispi;
1056 double * profile = NULL;
1057 const cpl_size i0 = 0;
1058 const double p0 = cpl_polynomial_get_coeff(disp, &i0);
1060 double xpos = (double)(1-hsize)-xtrunc;
1061 const double xmax = (double)(nself-hsize)+xtrunc;
1062 double xderiv, xextreme;
1063 cpl_error_code error = CPL_ERROR_NONE;
1065 cpl_size ulines = 0;
1067 cpl_ensure_code(
self != NULL, CPL_ERROR_NULL_INPUT);
1068 cpl_ensure_code(disp != NULL, CPL_ERROR_NULL_INPUT);
1069 cpl_ensure_code(lines != NULL, CPL_ERROR_NULL_INPUT);
1071 cpl_ensure_code(wslit > 0.0, CPL_ERROR_ILLEGAL_INPUT);
1072 cpl_ensure_code(wfwhm > 0.0, CPL_ERROR_ILLEGAL_INPUT);
1073 cpl_ensure_code(hsize >= 0, CPL_ERROR_ILLEGAL_INPUT);
1074 cpl_ensure_code(xtrunc > 0.0, CPL_ERROR_ILLEGAL_INPUT);
1075 cpl_ensure_code(nself > 2 * hsize, CPL_ERROR_ILLEGAL_INPUT);
1077 cpl_ensure_code(cpl_polynomial_get_dimension(disp) == 1,
1078 CPL_ERROR_ILLEGAL_INPUT);
1079 cpl_ensure_code(cpl_polynomial_get_degree(disp) > 0,
1080 CPL_ERROR_ILLEGAL_INPUT);
1083 wl = cpl_polynomial_eval_1d(disp, xpos, &xderiv);
1085 if (wl <= 0.0)
return
1086 cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_INPUT, __FILE__,
1087 __LINE__,
"Non-positive wavelength at x=%g: "
1088 "P(x)=%g, P'(x)=%g", xpos, wl, xderiv);
1090 if (xderiv <= 0.0)
return
1091 cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_INPUT, __FILE__,
1092 __LINE__,
"Non-increasing dispersion at "
1093 "x=%g: P'(x)=%g, P(x)=%g", xpos, xderiv, wl);
1096 iline = cpl_vector_find(xlines, wl);
1099 if (dxlines[iline] < wl) iline++;
1101 if (iline >= nlines)
return
1102 cpl_error_set_message_macro(cpl_func, CPL_ERROR_DATA_NOT_FOUND, __FILE__,
1103 __LINE__,
"The %d-line catalogue has only "
1104 "lines below P(%g)=%g > %g", nlines, xpos,
1105 wl, dxlines[nlines-1]);
1107 memset(dself, 0, nself *
sizeof(
double));
1109 dispi = cpl_polynomial_duplicate(disp);
1112 cpl_polynomial_derivative(dispi, 0);
1114 prestate = cpl_errorstate_get();
1116 if (cpl_polynomial_solve_1d(dispi, 0.5*(nlines+1), &xextreme, 1)) {
1117 cpl_errorstate_set(prestate);
1118 }
else if (xpos < xextreme && xextreme < xmax) {
1119 cpl_polynomial_delete(dispi);
1120 return cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
1121 __FILE__, __LINE__,
"Non-monotone "
1122 "dispersion at x=%g: P'(x)=0, "
1123 "P(x)=%g", xextreme,
1124 cpl_polynomial_eval_1d(disp, xextreme,
1129 const int npix = 1+(int)xtrunc;
1131 if (erftmp != NULL && cpl_vector_get_size(erftmp) == npix &&
1132 cpl_vector_get(erftmp, 0) > 0.0) {
1133 profile = cpl_vector_get_data(erftmp);
1136 const double yval = 0.5 / wslit;
1137 const double x0p = 0.5 * wslit + 0.5;
1138 const double x0n = -0.5 * wslit + 0.5;
1144 if (erftmp == NULL) {
1145 profile = (
double*)cpl_malloc(
sizeof(
double)*(size_t)npix);
1147 cpl_vector_set_size(erftmp, npix);
1148 profile = cpl_vector_get_data(erftmp);
1151 profile[0] = 2.0 * yval * x1diff;
1153 for (ipix = 1; ipix < npix; ipix++) {
1154 const double x1 = (double)ipix;
1155 const double x1p = x1 + 0.5 * wslit + 0.5;
1156 const double x1n = x1 - 0.5 * wslit + 0.5;
1157 const double x0diff = x1diff;
1162 profile[ipix] = yval * (x1diff - x0diff);
1168 cpl_polynomial_copy(dispi, disp);
1174 xpos -= (wl - dxlines[iline]) / xderiv;
1177 for (; !error && iline < nlines; iline++) {
1181 if (dylines[iline] <= 0.0)
continue;
1184 if (plinepix != NULL && plinepix[iline] > 0.0) xpos = plinepix[iline];
1186 if (xpos > xmax) xpos = xmax;
1189 error = cpl_polynomial_set_coeff(dispi, &i0, p0 - dxlines[iline]) ||
1190 cpl_polynomial_solve_1d(dispi, xpos, &xpos, 1);
1195 cpl_msg_debug(cpl_func,
"Stopping spectrum fill at line %d/%d "
1196 "at xpos=%g > xmax=%g",
1197 iline, nlines, xpos, xmax);
1198 cpl_errorstate_dump(prestate, CPL_FALSE,
1200 cpl_errorstate_set(prestate);
1204 if (linepix != NULL && ulines) (void)cpl_vector_fill(linepix, 0.0);
1205 (void)cpl_error_set_message_macro(cpl_func, cpl_error_get_code(),
1207 "Could not find pixel-position "
1208 "of line %d/%d at wavelength=%g."
1209 " xpos=%g, xmax=%g",
1210 iline, nlines, dxlines[iline],
1213 }
else if (dofast) {
1214 const double frac = fabs(xpos - floor(xpos));
1215 #ifdef IRPLIB_WAVECAL_FAST_FAST
1216 const double frac0 = 1.0 - frac;
1226 = (en1nw - en1pw) / (ep1pw - en1pw - ep1nw + en1nw);
1229 const double frac1 = 1.0 - frac0;
1230 const double yval0 = frac0 * dylines[iline];
1231 const double yval1 = frac1 * dylines[iline];
1232 const int npix = 1+(int)xtrunc;
1234 int i0n = hsize - 1 + floor(xpos);
1238 cpl_boolean didline = CPL_FALSE;
1242 if (plinepix != NULL) plinepix[iline] = xpos;
1245 (void)cpl_error_set_message_macro(cpl_func,
1246 CPL_ERROR_UNSPECIFIED,
1248 "Illegal split at x=%g: %g + "
1249 "%g = 1", xpos, frac0, frac1);
1250 #ifdef IRPLIB_WAVEVAL_DEBUG
1252 cpl_msg_warning(cpl_func,
"profile split at x=%g: %g + %g = 1",
1253 xpos, frac0, frac1);
1257 for (ipix = 0; ipix < npix; ipix++, i0n--, i0p++, i1n--, i1p++) {
1259 if (i0n >= 0 && i0n < nself) {
1260 dself[i0n] += yval0 * profile[ipix];
1263 if (i1n >= 0 && i1n < nself && ipix + 1 < npix) {
1264 dself[i1n] += yval1 * profile[ipix+1];
1268 if (ipix == 0)
continue;
1270 if (i0p >= 0 && i0p < nself) {
1271 dself[i0p] += yval0 * profile[ipix];
1274 if (i1p >= 0 && i1p < nself && ipix + 1 < npix) {
1275 dself[i1p] += yval1 * profile[ipix+1];
1280 if (didline) ulines++;
1283 const double yval = 0.5 * dylines[iline] / wslit;
1284 const int ifirst = IRPLIB_MAX((
int)(xpos-xtrunc+0.5), 1-hsize);
1285 const int ilast = IRPLIB_MIN((
int)(xpos+xtrunc), nself-hsize);
1287 const double x0 = (double)ifirst - xpos;
1288 const double x0p = x0 + 0.5*wslit - 0.5;
1289 const double x0n = x0 - 0.5*wslit - 0.5;
1295 if (plinepix != NULL) plinepix[iline] = xpos;
1297 if (ilast >= ifirst) ulines++;
1299 for (ipix = ifirst; ipix <= ilast; ipix++) {
1300 const double x1 = (double)ipix - xpos;
1301 const double x1p = x1 + 0.5*wslit + 0.5;
1302 const double x1n = x1 - 0.5*wslit + 0.5;
1303 const double x0diff = x1diff;
1308 dself[ipix+hsize-1] += yval * (x1diff - x0diff);
1314 cpl_polynomial_delete(dispi);
1315 if (erftmp == NULL) cpl_free(profile);
1317 cpl_ensure_code(!error, cpl_error_get_code());
1321 for (i = 0; i < nself; i++) {
1322 dself[i] = dself[i] > 0.0 ? log(1.0 + dself[i]) : 0.0;
1327 cpl_error_set_message_macro(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
1328 __FILE__, __LINE__,
"The %d-line "
1329 "catalogue has no lines in the range "
1330 "%g -> P(%g)=%g", nlines, wl, xmax,
1331 cpl_polynomial_eval_1d(disp, xmax, NULL));
1333 if (pulines != NULL) *pulines = ulines;
1335 return CPL_ERROR_NONE;
1350 return x * erf( x / (sigma * CPL_MATH_SQRT2))
1351 + 2.0 * sigma/CPL_MATH_SQRT2PI * exp(-0.5 * x * x / (sigma * sigma));
1365 static double irplib_gsl_correlation(
const gsl_vector *
self,
void * data)
1368 irplib_multimin * mindata = (irplib_multimin *)data;
1369 cpl_errorstate prestate = cpl_errorstate_get();
1370 int nobs, nmodel, ndiff;
1373 cpl_ensure(
self != NULL, CPL_ERROR_NULL_INPUT, GSL_NAN);
1374 cpl_ensure(data != NULL, CPL_ERROR_NULL_INPUT, GSL_NAN);
1376 cpl_ensure(mindata->filler != NULL, CPL_ERROR_NULL_INPUT, GSL_NAN);
1377 cpl_ensure(mindata->observed != NULL, CPL_ERROR_NULL_INPUT, GSL_NAN);
1378 cpl_ensure(mindata->spectrum != NULL, CPL_ERROR_NULL_INPUT, GSL_NAN);
1380 nobs = cpl_vector_get_size(mindata->observed);
1381 nmodel = cpl_vector_get_size(mindata->spectrum);
1382 ndiff = nmodel - nobs;
1384 cpl_ensure((ndiff & 1) == 0, CPL_ERROR_ILLEGAL_INPUT, GSL_NAN);
1386 cpl_ensure(cpl_vector_get_size(mindata->vxc) == 1 + ndiff,
1387 CPL_ERROR_ILLEGAL_INPUT, GSL_NAN);
1391 for (i=0; i < (cpl_size)self->size; i++) {
1392 const double value = gsl_vector_get(
self, (
size_t)i);
1393 cpl_polynomial_set_coeff(mindata->disp1d, &i, value);
1398 cpl_ensure_code(!cpl_polynomial_shift_1d(mindata->disp1d, 0, -ndiff),
1399 cpl_error_get_code());
1401 if (mindata->filler(mindata->spectrum, mindata->disp1d,
1403 || !cpl_errorstate_is_equal(prestate)) {
1408 (void)cpl_vector_fill(mindata->vxc, -1.0);
1410 mindata->maxxc = ndiff;
1412 if (!cpl_errorstate_is_equal(prestate)) {
1413 cpl_msg_debug(cpl_func,
"Spectrum fill failed:");
1414 cpl_errorstate_dump(prestate, CPL_FALSE,
1416 cpl_errorstate_set(prestate);
1420 mindata->maxxc = cpl_vector_correlate(mindata->vxc,
1425 #ifdef IRPLIB_SPC_DUMP
1427 irplib_polynomial_dump_corr_step(mindata->disp1d, mindata->vxc,
1431 mindata->xc = cpl_vector_get(mindata->vxc, ndiff);
1433 if (mindata->maxxc != ndiff &&
1434 cpl_vector_get(mindata->vxc, mindata->maxxc) > mindata->mxc) {
1435 const irplib_base_spectrum_model * arclamp
1436 = (
const irplib_base_spectrum_model *)mindata->param;
1438 if (mindata->mdisp == NULL) {
1439 mindata->mdisp = cpl_polynomial_duplicate(mindata->disp1d);
1441 cpl_polynomial_copy(mindata->mdisp, mindata->disp1d);
1443 mindata->mxc = cpl_vector_get(mindata->vxc, mindata->maxxc);
1444 mindata->ishift = mindata->maxxc;
1445 cpl_msg_debug(cpl_func,
"Local maximum: %g(%d) > %g(%d) (cost=%u:%u. "
1446 "lines=%u)", mindata->mxc, mindata->maxxc, mindata->xc,
1447 ndiff, (
unsigned)arclamp->cost, (
unsigned)arclamp->xcost,
1448 (
unsigned)arclamp->ulines);
1451 return -mindata->xc;
1483 const cpl_vector * obs,
1486 irplib_base_spectrum_model* model,
1487 cpl_error_code (* filler)
1489 const cpl_polynomial *,
1490 irplib_base_spectrum_model *),
1503 cpl_errorstate prestate = cpl_errorstate_get();
1504 cpl_polynomial * start;
1505 cpl_polynomial * cand;
1506 cpl_polynomial * backup;
1507 cpl_error_code error = CPL_ERROR_NONE;
1509 cpl_bivector * xtshift = cpl_bivector_new(nmaxima ? nmaxima : 1);
1510 const cpl_vector * xtshiftx = cpl_bivector_get_x_const(xtshift);
1511 const cpl_vector * xtshifty = cpl_bivector_get_y_const(xtshift);
1518 cpl_ensure_code(
self != NULL, CPL_ERROR_NULL_INPUT);
1519 cpl_ensure_code(obs != NULL, CPL_ERROR_NULL_INPUT);
1520 cpl_ensure_code(model != NULL, CPL_ERROR_NULL_INPUT);
1521 cpl_ensure_code(filler != NULL, CPL_ERROR_NULL_INPUT);
1522 cpl_ensure_code(pxc != NULL, CPL_ERROR_NULL_INPUT);
1524 cpl_ensure_code(cpl_polynomial_get_dimension(
self) == 1,
1525 CPL_ERROR_ILLEGAL_INPUT);
1527 cpl_ensure_code(cpl_polynomial_get_degree(
self) > 0,
1528 CPL_ERROR_ILLEGAL_INPUT);
1530 cpl_ensure_code(maxdeg >= 0, CPL_ERROR_ILLEGAL_INPUT);
1531 cpl_ensure_code(pixtol > 0.0, CPL_ERROR_ILLEGAL_INPUT);
1532 cpl_ensure_code(pixstep > 0.0, CPL_ERROR_ILLEGAL_INPUT);
1533 cpl_ensure_code(hsize >= 0, CPL_ERROR_ILLEGAL_INPUT);
1534 cpl_ensure_code(maxite >= 0, CPL_ERROR_ILLEGAL_INPUT);
1535 cpl_ensure_code(nmaxima >= 0, CPL_ERROR_ILLEGAL_INPUT);
1536 cpl_ensure_code(maxfail > 0, CPL_ERROR_ILLEGAL_INPUT);
1537 cpl_ensure_code(maxcont > 0, CPL_ERROR_ILLEGAL_INPUT);
1538 cpl_ensure_code(linelim >= 0, CPL_ERROR_ILLEGAL_INPUT);
1542 cpl_ensure_code(doplot == CPL_TRUE || doplot == CPL_FALSE,
1543 CPL_ERROR_ILLEGAL_INPUT);
1544 return cpl_error_set_message(cpl_func, CPL_ERROR_UNSUPPORTED_MODE,
1545 "GSL is not available");
1550 hsize, doplot, &xc)) {
1551 cpl_bivector_delete(xtshift);
1552 return cpl_error_set_where(cpl_func);
1555 if (model->ulines > (cpl_size)linelim) {
1557 const double xxc = cpl_vector_get(xtshiftx, 0);
1558 const double xc0 = cpl_vector_get(xtshifty, 0);
1560 cpl_msg_warning(cpl_func,
"Doing only shift=%g pixels with lines=%u > "
1561 "%d and XC=%g", xxc, (
unsigned)model->ulines, linelim,
1564 cpl_polynomial_shift_1d(
self, 0, xxc);
1568 cpl_bivector_delete(xtshift);
1570 return CPL_ERROR_NONE;
1573 start = cpl_polynomial_duplicate(
self);
1574 cand = cpl_polynomial_new(1);
1575 backup = cpl_polynomial_new(1);
1578 nshift = cpl_bivector_get_size(xtshift);
1579 if (nmaxima == 0 || nmaxima > nshift) nmaxima = nshift;
1581 cpl_msg_info(cpl_func,
"Optimizing %d/%d local shift-maxima "
1582 "(no-shift xc=%g. linelim=%d)", nmaxima, nshift, xc, linelim);
1583 if (cpl_msg_get_level() <= CPL_MSG_DEBUG)
1584 cpl_bivector_dump(xtshift, stdout);
1586 for (imaxima = 0; imaxima < nmaxima; imaxima++) {
1588 const double xxc = cpl_vector_get(xtshiftx, imaxima);
1589 double xtpixstep = pixstep;
1590 double xtpixtol = pixtol;
1592 cpl_boolean ok = CPL_FALSE;
1596 cpl_polynomial_copy(cand, start);
1597 cpl_polynomial_shift_1d(cand, 0, xxc);
1598 cpl_polynomial_copy(backup, cand);
1601 for (nfail = 0; nfail < maxfail; nfail++, xtpixtol *= 2.0,
1603 int restart = maxcont;
1608 cpl_errorstate_dump(prestate, CPL_FALSE,
1610 cpl_errorstate_set(prestate);
1612 error = irplib_polynomial_find_1d_from_correlation_
1613 (cand, maxdeg, obs, model,
1614 filler, xtpixtol, xtpixstep, 2,
1615 maxite, &xtxc, &redo);
1616 if (redo && !error) error = CPL_ERROR_CONTINUE;
1617 }
while (((!error && redo) || error == CPL_ERROR_CONTINUE)
1620 if (!error && !redo) {
1621 cpl_msg_debug(cpl_func,
"XC(imax=%d/%d:xtpixtol=%g): %g "
1622 "(cost=%u:%u)", 1+imaxima, nmaxima, xtpixtol,
1623 xtxc, (
unsigned)model->cost,
1624 (
unsigned)model->xcost);
1627 cpl_msg_warning(cpl_func,
"Increasing xtpixtol from %g (%g, imax="
1628 "%d/%d)", xtpixtol, xtpixstep, 1+imaxima, nmaxima);
1629 if (model->ulines > (cpl_size)linelim) {
1630 cpl_msg_warning(cpl_func,
"Stopping search-refinement via "
1631 "catalogue with %u lines > %d",
1632 (
unsigned)model->ulines, linelim);
1635 cpl_polynomial_copy(cand, start);
1639 for (; !error && xtpixtol > 0.0; xtpixtol *= 0.25, xtpixstep *= 0.5) {
1640 int restart = maxcont;
1643 cpl_polynomial_copy(backup, cand);
1646 cpl_errorstate_dump(prestate, CPL_FALSE,
1648 cpl_errorstate_set(prestate);
1650 error = irplib_polynomial_find_1d_from_correlation_
1651 (cand, maxdeg, obs, model, filler,
1652 xtpixtol, xtpixstep, 2, maxite, &xtxc, &redo);
1653 if (redo && !error) error = CPL_ERROR_CONTINUE;
1654 }
while (((!error && redo) || error == CPL_ERROR_CONTINUE)
1659 cpl_msg_debug(cpl_func,
"XC(imax=%d/%d:xtpixtol=%g): %g (cost=%u:%u"
1660 ". ulines=%u)", 1+imaxima, nmaxima, xtpixtol, xtxc,
1661 (
unsigned)model->cost, (
unsigned)model->xcost,
1662 (
unsigned)model->ulines);
1663 if (model->ulines > (cpl_size)linelim) {
1664 cpl_msg_info(cpl_func,
"Stopping search-refinement via "
1665 "catalogue with %u lines > %u",
1666 (
unsigned)model->ulines, linelim);
1673 cpl_errorstate_dump(prestate, CPL_FALSE,
1675 cpl_errorstate_set(prestate);
1676 cpl_polynomial_copy(cand, backup);
1678 if (ok && xtxc > xc) {
1680 cpl_polynomial_copy(
self, cand);
1683 cpl_msg_info(cpl_func,
"XC(imax=%d/%d): %g -> %g (initial-shift=%g. "
1684 "cost=%u:%u. lines=%u)", 1+imaxima, nmaxima,
1685 cpl_vector_get(xtshifty, imaxima), xtxc,
1686 cpl_vector_get(xtshiftx, imaxima),
1687 (
unsigned)model->cost, (
unsigned)model->xcost,
1688 (
unsigned)model->ulines);
1690 cpl_msg_info(cpl_func,
"xc(imax=%d/%d): %g -> %g (initial-shift=%g. "
1691 "cost=%u:%u. lines=%u)", 1+imaxima, nmaxima,
1692 cpl_vector_get(xtshifty, imaxima), xtxc,
1693 cpl_vector_get(xtshiftx, imaxima),
1694 (
unsigned)model->cost, (
unsigned)model->xcost,
1695 (
unsigned)model->ulines);
1699 cpl_polynomial_delete(start);
1700 cpl_polynomial_delete(backup);
1701 cpl_polynomial_delete(cand);
1705 const double xxc = cpl_vector_get(xtshiftx, 0);
1706 const double xc0 = cpl_vector_get(xtshifty, 0);
1708 error = cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
1709 "Could not improve XC=%g over %d "
1710 "local shift-maxima, best at shift %g",
1713 cpl_msg_info(cpl_func,
"Maximal XC=%g (up from %g, with initial pixel-"
1714 "shift of %g) at %d/%d local shift-maximi", xc,
1715 cpl_vector_get(xtshifty, imaximum),
1716 cpl_vector_get(xtshiftx, imaximum),
1717 1+imaximum, nmaxima);
1726 cpl_bivector_delete(xtshift);
cpl_error_code irplib_vector_fill_line_spectrum_model(cpl_vector *self, cpl_vector *linepix, cpl_vector *erftmp, const cpl_polynomial *disp, const cpl_bivector *lines, double wslit, double wfwhm, double xtrunc, int hsize, cpl_boolean dofast, cpl_boolean dolog, cpl_size *pulines)
Generate a 1D spectrum from (arc) lines and a dispersion relation.
cpl_error_code irplib_vector_fill_line_spectrum(cpl_vector *self, const cpl_polynomial *disp, irplib_base_spectrum_model *lsslamp)
Generate a 1D spectrum from a model and a dispersion relation.
cpl_error_code irplib_bivector_find_shift_from_correlation(cpl_bivector *self, const cpl_polynomial *disp, const cpl_vector *obs, irplib_base_spectrum_model *model, cpl_error_code(*filler)(cpl_vector *, const cpl_polynomial *, irplib_base_spectrum_model *), int hsize, cpl_boolean doplot, double *pxc)
Find shift(s) that maximizes (locally) the cross-correlation.
cpl_error_code irplib_vector_fill_logline_spectrum_fast(cpl_vector *self, const cpl_polynomial *disp, irplib_base_spectrum_model *lsslamp)
Generate a 1D spectrum from a model and a dispersion relation.
cpl_error_code irplib_polynomial_shift_1d_from_correlation(cpl_polynomial *self, const cpl_vector *obs, irplib_base_spectrum_model *model, cpl_error_code(*filler)(cpl_vector *, const cpl_polynomial *, irplib_base_spectrum_model *), int hsize, cpl_boolean doplot, double *pxc)
Shift self by the amount that maximizes the cross-correlation.
cpl_error_code irplib_polynomial_find_1d_from_correlation(cpl_polynomial *self, int maxdeg, const cpl_vector *obs, irplib_base_spectrum_model *model, cpl_error_code(*filler)(cpl_vector *, const cpl_polynomial *, irplib_base_spectrum_model *), double pixtol, double pixstep, int hsize, int maxite, double *pxc)
Modify self by maximizing the cross-correlation.
void irplib_errorstate_dump_debug(unsigned self, unsigned first, unsigned last)
Dump a single CPL error at the CPL debug level.
cpl_error_code irplib_vector_fill_logline_spectrum(cpl_vector *self, const cpl_polynomial *disp, irplib_base_spectrum_model *lsslamp)
Generate a 1D spectrum from a model and a dispersion relation.
cpl_error_code irplib_polynomial_fit_2d_dispersion(cpl_polynomial *self, const cpl_image *imgwave, int fitdeg, double *presid)
Fit a 2D-dispersion from an image of wavelengths.
cpl_error_code irplib_vector_fill_line_spectrum_fast(cpl_vector *self, const cpl_polynomial *disp, irplib_base_spectrum_model *lsslamp)
Generate a 1D spectrum from a model and a dispersion relation.
cpl_error_code irplib_plot_spectrum_and_model(const cpl_vector *self, const cpl_polynomial *disp1d, irplib_base_spectrum_model *model, cpl_error_code(*filler)(cpl_vector *, const cpl_polynomial *, irplib_base_spectrum_model *))
Plot a 1D spectrum and one from a model.
cpl_error_code irplib_polynomial_find_1d_from_correlation_all(cpl_polynomial *self, int maxdeg, const cpl_vector *obs, int nmaxima, int linelim, irplib_base_spectrum_model *model, cpl_error_code(*filler)(cpl_vector *, const cpl_polynomial *, irplib_base_spectrum_model *), double pixtol, double pixstep, int hsize, int maxite, int maxfail, int maxcont, cpl_boolean doplot, double *pxc)
Modify self by maximizing the cross-correlation across all maxima.
int irplib_bivector_count_positive(const cpl_bivector *self, double x_min, double x_max)
Count the positive Y-entries in a given X-range.
double irplib_erf_antideriv(double x, double sigma)
The antiderivative of erx(x/sigma/sqrt(2)) with respect to x.