36#include "irplib_strehl.h"
37#include "irplib_utils.h"
53#ifndef IRPLIB_STREHL_RAD_CENTRAL
54#define IRPLIB_STREHL_RAD_CENTRAL 5
57#ifndef IRPLIB_STREHL_DETECT_LEVEL
58#define IRPLIB_STREHL_DETECT_LEVEL 5.0
61#define IRPLIB_DISK_BG_MIN_PIX_NB 30
62#define IRPLIB_DISK_BG_REJ_LOW 0.1
63#define IRPLIB_DISK_BG_REJ_HIGH 0.1
66#define IRPLIB_MIN CPL_MIN
68#define IRPLIB_MIN(A,B) (((A) < (B)) ? (A) : (B))
72#define IRPLIB_MAX CPL_MAX
74#define IRPLIB_MAX(A,B) (((A) > (B)) ? (A) : (B))
81static cpl_image * irplib_strehl_generate_otf(
double,
double,
double,
double,
83static double PSF_H1(
double,
double,
double);
84static double PSF_H2(
double,
double);
85static double PSF_G(
double,
double);
86static double PSF_sinc_norm(
double);
87static double PSF_TelOTF(
double,
double);
89#ifndef IRPLIB_NO_FIT_GAUSSIAN
90#ifdef IRPLIB_STREHL_USE_CPL_IMAGE_FIT_GAUSSIAN
91static double irplib_gaussian_2d(
double,
double,
double,
double,
double);
94#if defined CPL_VERSION_CODE && CPL_VERSION_CODE >= CPL_VERSION(6, 9, 1)
95#define irplib_gaussian_eval_2d cpl_gaussian_eval_2d
97static double irplib_gaussian_eval_2d(
const cpl_array *,
double,
double);
100static uint32_t irplib_roundup_power2(uint32_t v) CPL_ATTR_CONST;
103cpl_error_code irplib_gaussian_maxpos(
const cpl_image *,
112irplib_closeset_aperture(
const cpl_apertures * self,
113 const double x,
const double y,
int * ind);
152cpl_error_code irplib_strehl_compute(
const cpl_image * im,
176 double star_radius, max_radius;
179 const double window_size = (double)(IRPLIB_STREHL_RAD_CENTRAL);
182 const double strehl_error_coefficient = CPL_MATH_PI * 0.007 / 0.0271;
186#ifndef IRPLIB_NO_FIT_GAUSSIAN
187 double xposfit = 0.0, yposfit = 0.0, peak = 0.0;
190 cpl_errorstate prestate = cpl_errorstate_get();
193 cpl_ensure_code(window_size > 0.0, CPL_ERROR_ILLEGAL_INPUT);
196 cpl_ensure_code(im != NULL, CPL_ERROR_NULL_INPUT);
197 cpl_ensure_code(strehl != NULL, CPL_ERROR_NULL_INPUT);
198 cpl_ensure_code(strehl_err != NULL, CPL_ERROR_NULL_INPUT);
199 cpl_ensure_code(star_bg != NULL, CPL_ERROR_NULL_INPUT);
200 cpl_ensure_code(star_peak != NULL, CPL_ERROR_NULL_INPUT);
201 cpl_ensure_code(star_flux != NULL, CPL_ERROR_NULL_INPUT);
202 cpl_ensure_code(psf_peak != NULL, CPL_ERROR_NULL_INPUT);
203 cpl_ensure_code(psf_flux != NULL, CPL_ERROR_NULL_INPUT);
205 cpl_ensure_code(pscale > 0.0, CPL_ERROR_ILLEGAL_INPUT);
207 cpl_ensure_code(r1 > 0.0, CPL_ERROR_ILLEGAL_INPUT);
208 cpl_ensure_code(r2 > 0.0, CPL_ERROR_ILLEGAL_INPUT);
209 cpl_ensure_code(r3 > r2, CPL_ERROR_ILLEGAL_INPUT);
215 psf = irplib_strehl_generate_psf(m1, m2, lam, dlam, pscale, size);
217 return cpl_error_set_where(cpl_func);
221 *psf_peak = cpl_image_get_max(psf);
222 cpl_image_delete(psf);
224 assert( *psf_peak > 0.0);
227#ifndef IRPLIB_NO_FIT_GAUSSIAN
228 code = irplib_gaussian_maxpos(im, IRPLIB_STREHL_DETECT_LEVEL, xpos, ypos,
229 &xposfit, &yposfit, &peak);
231 cpl_errorstate_set(prestate);
239 *star_bg = irplib_strehl_ring_background(im, xpos, ypos,
240 r2/pscale, r3/pscale,
241 IRPLIB_BG_METHOD_AVER_REJ);
242 if (!cpl_errorstate_is_equal(prestate)) {
243 return cpl_error_set_where(cpl_func);
247 star_radius = r1/pscale;
250 *star_flux = irplib_strehl_disk_flux(im, xpos, ypos, star_radius, *star_bg);
252 if (*star_flux <= 0.0) {
253 return cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_OUTPUT,
254 "Non-positive star flux=%g (Star "
255 "background=%g)", *star_flux, *star_bg);
259 max_radius = window_size < star_radius ? window_size : star_radius;
260 cpl_ensure_code(!irplib_strehl_disk_max(im, xpos, ypos, max_radius,
261 star_peak), cpl_error_get_code());
262 *star_peak -= *star_bg;
264 if (*star_flux <= 0.0) {
265 return cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_OUTPUT,
266 "Non-positive star peak=%g (Star "
267 "background=%g, Star flux=%g)",
268 *star_flux, *star_bg, *star_flux);
273 *strehl = (*star_peak * *psf_flux ) / ( *star_flux * *psf_peak);
275#ifndef IRPLIB_NO_FIT_GAUSSIAN
276 if (code == CPL_ERROR_NONE && peak > *star_peak && *star_peak > 0.0 &&
277 *strehl * peak / *star_peak <= 1.0) {
278 cpl_msg_debug(cpl_func,
"Increasing Strehl from %g: %g (%g)",
279 *strehl, *strehl * peak / *star_peak,
281 *strehl *= peak / *star_peak;
292 while (cpl_flux_get_noise_ring(im, ring, noise_box_sz, noise_nsamples,
293 bg_noise, NULL) && --ring_tries > 0);
294 if (ring_tries > 0) {
295 cpl_errorstate_set(prestate);
297 return cpl_error_set_where(cpl_func);
300 *strehl_err = strehl_error_coefficient * (*bg_noise) * pscale *
301 star_radius * star_radius / *star_flux;
304 cpl_msg_warning(cpl_func,
"Extreme Strehl-ratio=%g (strehl-error=%g, "
305 "star_peak=%g, star_flux=%g, psf_peak=%g, psf_flux=%g)",
306 *strehl, *strehl_err, *star_peak, *star_flux, *psf_peak,
311 return *strehl_err >= 0.0
313 : cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_OUTPUT,
314 "Negative strehl-error=%g (Strehl-ratio=%g, "
315 "star_peak=%g, star_flux=%g, psf_peak=%g, "
316 "psf_flux=%g", *strehl_err, *strehl,
317 *star_peak, *star_flux, *psf_peak, *psf_flux);
334double irplib_strehl_disk_flux(
const cpl_image * im,
340 const int nx = cpl_image_get_size_x(im);
341 const int ny = cpl_image_get_size_y(im);
343 const int lx = (int)(xpos - rad);
344 const int ly = (int)(ypos - rad);
346 const int ux = (int)(xpos + rad) + 1;
347 const int uy = (int)(ypos + rad) + 1;
349 const double sqr = rad * rad;
355 cpl_ensure(im != NULL, CPL_ERROR_NULL_INPUT, 0.0);
356 cpl_ensure(rad > 0.0, CPL_ERROR_ILLEGAL_INPUT, 0.0);
358 for (j = IRPLIB_MAX(ly, 0); j < IRPLIB_MIN(uy, ny); j++) {
359 const double yj = (double)j - ypos;
360 for (i = IRPLIB_MAX(lx, 0); i < IRPLIB_MIN(ux, nx); i++) {
361 const double xi = (double)i - xpos;
362 const double dist = yj * yj + xi * xi;
365 const double value = cpl_image_get(im, i+1, j+1, &isbad);
392double irplib_strehl_ring_background(
const cpl_image * im,
397 irplib_strehl_bg_method mode)
399 const int nx = cpl_image_get_size_x(im);
400 const int ny = cpl_image_get_size_y(im);
402 const int lx = (int)(xpos - rad_ext);
403 const int ly = (int)(ypos - rad_ext);
405 const int ux = (int)(xpos + rad_ext) + 1;
406 const int uy = (int)(ypos + rad_ext) + 1;
408 const double sqr_int = rad_int * rad_int;
409 const double sqr_ext = rad_ext * rad_ext;
410 cpl_vector * pix_arr;
415 cpl_ensure(im != NULL, CPL_ERROR_NULL_INPUT, 0.0);
416 cpl_ensure(rad_int > 0.0, CPL_ERROR_ILLEGAL_INPUT, 0.0);
417 cpl_ensure(rad_ext > rad_int, CPL_ERROR_ILLEGAL_INPUT, 0.0);
419 cpl_ensure(mode == IRPLIB_BG_METHOD_AVER_REJ ||
420 mode == IRPLIB_BG_METHOD_MEDIAN,
421 CPL_ERROR_UNSUPPORTED_MODE, 0.0);
423 mpix = (int)((2.0 * rad_ext + 1.0) * (2.0 * rad_ext + 1.0));
426 pix_arr = cpl_vector_new(mpix);
431 for (j = IRPLIB_MAX(ly, 0); j < IRPLIB_MIN(uy, ny); j++) {
432 const double yj = (double)j - ypos;
433 for (i = IRPLIB_MAX(lx, 0); i < IRPLIB_MIN(ux, nx); i++) {
434 const double xi = (double)i - xpos;
435 const double dist = yj * yj + xi * xi;
436 if (sqr_int <= dist && dist <= sqr_ext) {
438 const double value = cpl_image_get(im, i+1, j+1, &isbad);
441 cpl_vector_set(pix_arr, npix, value);
448 assert(npix <= mpix);
450 if (npix < IRPLIB_DISK_BG_MIN_PIX_NB) {
451 cpl_vector_delete(pix_arr);
452 (void)cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
"Need "
453 "at least %d (not %d <= %d) samples to "
454 "compute noise", IRPLIB_DISK_BG_MIN_PIX_NB,
462 pix_arr = cpl_vector_wrap(npix, (
double*)cpl_vector_unwrap(pix_arr));
464 if (mode == IRPLIB_BG_METHOD_AVER_REJ) {
465 const int low_ind = (int)((
double)npix * IRPLIB_DISK_BG_REJ_LOW);
466 const int high_ind = (int)((
double)npix
467 * (1.0 - IRPLIB_DISK_BG_REJ_HIGH));
470 cpl_vector_sort(pix_arr, CPL_SORT_ASCENDING);
472 for (i=low_ind; i<high_ind; i++) {
473 flux += cpl_vector_get(pix_arr, i);
475 if (high_ind - low_ind > 1) flux /= (double)(high_ind - low_ind);
477 flux = cpl_vector_get_median(pix_arr);
480 cpl_vector_delete(pix_arr);
506cpl_image * irplib_strehl_generate_psf(
double m1,
513 cpl_image * otf_image = irplib_strehl_generate_otf(m1, m2, lam, dlam,
516 if (otf_image == NULL ||
527 cpl_image_fft(otf_image, NULL, CPL_FFT_UNNORMALIZED) ||
530 cpl_image_abs(otf_image) ||
533 cpl_image_normalise(otf_image, CPL_NORM_FLUX)) {
535 (void)cpl_error_set_where(cpl_func);
536 cpl_image_delete(otf_image);
562static cpl_image * irplib_strehl_generate_otf(
double m1,
571 const double obs_ratio = m1 != 0.0 ? m2 / m1 : 0.0;
573 const double rpscale = pscale * CPL_MATH_2PI / (double)(360 * 60 * 60);
575 const double f_max = m1 * rpscale * (double)size;
578 const int pix0 = size / 2;
582 cpl_ensure(m2 > 0.0, CPL_ERROR_ILLEGAL_INPUT, NULL);
583 cpl_ensure(m1 > m2, CPL_ERROR_ILLEGAL_INPUT, NULL);
584 cpl_ensure(dlam > 0.0, CPL_ERROR_ILLEGAL_INPUT, NULL);
585 cpl_ensure(pscale > 0.0, CPL_ERROR_ILLEGAL_INPUT, NULL);
586 cpl_ensure(size > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
588 cpl_ensure(size % 2 == 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
591 cpl_ensure(2.0 * lam > dlam, CPL_ERROR_ILLEGAL_INPUT, NULL);
598 otf_data = (
double*)cpl_malloc(size * size *
sizeof(*otf_data));
605 for (j = 0; j <= pix0; j++) {
606 double sinc_y_9 = 0.0;
607 for (i = 0; i <= j; i++) {
608 if (i == 0 && j == 0) {
609 otf_data[size * pix0 + pix0] = 1.0;
611 const double x = (double)i;
612 const double y = (double)j;
613 const double sqdist = x * x + y * y;
614 double f_lambda = 0.0;
615 double sinc_xy_9 = 0.0;
623 for (k = 4; k >= -4; k--) {
625 const double lambda = lam - dlam * (double)k / 8.0;
629 if (sqdist * lambda * lambda >= f_max * f_max)
break;
632 f_lambda = sqrt(sqdist) / f_max;
635 sinc_xy_9 = sinc_y_9 =
636 PSF_sinc_norm(y / (
double)size) / 9.0;
638 sinc_xy_9 = sinc_y_9 *
639 PSF_sinc_norm(x / (
double)size);
643 otfxy += PSF_TelOTF(f_lambda * lambda, obs_ratio);
649 otf_data[size * (pix0 - j) + pix0 - i] = otfxy;
650 otf_data[size * (pix0 - i) + pix0 - j] = otfxy;
652 otf_data[size * (pix0 - j) + pix0 + i] = otfxy;
653 otf_data[size * (pix0 + i) + pix0 - j] = otfxy;
655 otf_data[size * (pix0 + j) + pix0 - i] = otfxy;
656 otf_data[size * (pix0 - i) + pix0 + j] = otfxy;
657 otf_data[size * (pix0 + j) + pix0 + i] = otfxy;
658 otf_data[size * (pix0 + i) + pix0 + j] = otfxy;
665 return cpl_image_wrap_double(size, size, otf_data);
676 const double e = fabs(1.0-v) > 0.0 ? -1.0 : 1.0;
678 return((v*v/CPL_MATH_PI)*acos((f/v)*(1.0+e*(1.0-u*u)/(4.0*f*f))));
684static double PSF_H2(
double f,
687 const double tmp1 = (2.0 * f) / (1.0 + u);
688 const double tmp2 = (1.0 - u) / (2.0 * f);
690 return -1.0 * (f/CPL_MATH_PI) * (1.0+u)
691 * sqrt((1.0-tmp1*tmp1)*(1.0-tmp2*tmp2));
697static double PSF_G(
double f,
700 if (f <= (1.0-u)/2.0)
return(u*u);
701 if (f >= (1.0+u)/2.0)
return(0.0);
702 else return(PSF_H1(f,u,1.0) + PSF_H1(f,u,u) + PSF_H2(f,u));
714static double PSF_sinc_norm(
double x)
718 return x != 0.0 ? sin(x * CPL_MATH_PI) / (x * CPL_MATH_PI) : 1.0;
724static double PSF_TelOTF(
double f,
727 return((PSF_G(f,1.0)+u*u*PSF_G(f/u,1.0)-2.0*PSF_G(f,u))/(1.0-u*u));
742cpl_error_code irplib_strehl_disk_max(
const cpl_image * self,
749 const int nx = cpl_image_get_size_x(self);
750 const int ny = cpl_image_get_size_y(self);
752 const int lx = (int)(xpos - radius);
753 const int ly = (int)(ypos - radius);
755 const int ux = (int)(xpos + radius) + 1;
756 const int uy = (int)(ypos + radius) + 1;
758 const double sqr = radius * radius;
759 cpl_boolean first = CPL_TRUE;
764 cpl_ensure_code(self != NULL, CPL_ERROR_NULL_INPUT);
765 cpl_ensure_code(ppeak != NULL, CPL_ERROR_NULL_INPUT);
766 cpl_ensure_code(radius > 0.0, CPL_ERROR_ILLEGAL_INPUT);
769 for (j = IRPLIB_MAX(ly, 0); j < IRPLIB_MIN(uy, ny); j++) {
770 const double yj = (double)j - ypos;
771 for (i = IRPLIB_MAX(lx, 0); i < IRPLIB_MIN(ux, nx); i++) {
772 const double xi = (double)i - xpos;
773 const double dist = yj * yj + xi * xi;
776 const double value = cpl_image_get(self, i+1, j+1, &isbad);
779 (first || value > *ppeak)) {
788 ? cpl_error_set(cpl_func, CPL_ERROR_DATA_NOT_FOUND)
792#ifndef IRPLIB_NO_FIT_GAUSSIAN
793#ifdef IRPLIB_STREHL_USE_CPL_IMAGE_FIT_GAUSSIAN
811static double irplib_gaussian_2d(
double x,
819 return norm / (sig_x * sig_y * CPL_MATH_2PI *
820 exp(x * x / (2.0 * sig_x * sig_x) +
821 y * y / (2.0 * sig_y * sig_y)));
825#if defined CPL_VERSION_CODE && CPL_VERSION_CODE >= CPL_VERSION(6, 9, 1)
848double irplib_gaussian_eval_2d(
const cpl_array * self,
double x,
double y)
850 cpl_errorstate prestate = cpl_errorstate_get();
851 const double B = cpl_array_get_double(self, 0, NULL);
852 const double A = cpl_array_get_double(self, 1, NULL);
853 const double R = cpl_array_get_double(self, 2, NULL);
854 const double M_x = cpl_array_get_double(self, 3, NULL);
855 const double M_y = cpl_array_get_double(self, 4, NULL);
856 const double S_x = cpl_array_get_double(self, 5, NULL);
857 const double S_y = cpl_array_get_double(self, 6, NULL);
861 if (!cpl_errorstate_is_equal(prestate)) {
862 (void)cpl_error_set_where(cpl_func);
863 }
else if (cpl_array_get_size(self) != 7) {
864 (void)cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_INPUT);
865 }
else if (fabs(R) < 1.0 && S_x != 0.0 && S_y != 0.0) {
866 const double x_n = (x - M_x) / S_x;
867 const double y_n = (y - M_y) / S_y;
869 value = B + A / (CPL_MATH_2PI * S_x * S_y * sqrt(1 - R * R)) *
870 exp(-0.5 / (1 - R * R) * ( x_n * x_n + y_n * y_n
871 - 2.0 * R * x_n * y_n));
872 }
else if (fabs(R) > 1.0) {
873 (void)cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_OUTPUT,
874 "fabs(R=%g) > 1", R);
876 (void)cpl_error_set_message(cpl_func, CPL_ERROR_DIVISION_BY_ZERO,
877 "R=%g. Sigma=(%g, %g)", R, S_x, S_y);
892static uint32_t irplib_roundup_power2(uint32_t v)
905irplib_closeset_aperture(
const cpl_apertures * self,
906 const double x,
const double y,
int * ind){
908 double dist = INFINITY;
910 const int nsize = cpl_apertures_get_size(self);
913 cpl_ensure_code(nsize > 0, cpl_error_get_code());
914 cpl_ensure_code(ind, CPL_ERROR_NULL_INPUT);
918 for(i = 1; i <= nsize; ++i){
919 const double this_x = cpl_apertures_get_centroid_x(self, i);
920 const double this_y = cpl_apertures_get_centroid_y(self, i);
921 const double this_d_sq = pow(this_x - x, 2.0) + pow(this_y - y, 2.0);
923 if(this_d_sq < dist){
929 return CPL_ERROR_NONE;
945cpl_error_code irplib_gaussian_maxpos(
const cpl_image * self,
954 const cpl_size nx = cpl_image_get_size_x(self);
955 const cpl_size ny = cpl_image_get_size_y(self);
959 const double median = cpl_image_get_median_dev(self, &med_dist);
960 cpl_mask * selection;
961 cpl_size nlabels = 0;
962 cpl_image * labels = NULL;
963 cpl_apertures * aperts;
967 cpl_size xposmax, yposmax;
968 double xposcen, yposcen;
969 double valmax, valfit = -1.0;
970#ifdef IRPLIB_STREHL_USE_CPL_IMAGE_FIT_GAUSSIAN
971 double norm, xcen, ycen, sig_x, sig_y, fwhm_x, fwhm_y;
973 cpl_array * gauss_parameters = NULL;
974 cpl_errorstate prestate = cpl_errorstate_get();
975 cpl_error_code code = CPL_ERROR_NONE;
978 cpl_ensure_code( sigma > 0.0, CPL_ERROR_ILLEGAL_INPUT);
980 selection = cpl_mask_new(nx, ny);
982 for (; iretry > 0 && nlabels == 0; iretry--, sigma *= 0.5) {
985 const double threshold = median + sigma * med_dist;
989 code = cpl_mask_threshold_image(selection, self, threshold, DBL_MAX,
995 cpl_image_delete(labels);
996 labels = cpl_image_labelise_mask_create(selection, &nlabels);
1000 cpl_mask_delete(selection);
1003 cpl_image_delete(labels);
1004 return cpl_error_set_where(cpl_func);
1005 }
else if (nlabels == 0) {
1006 cpl_image_delete(labels);
1007 return cpl_error_set(cpl_func, CPL_ERROR_DATA_NOT_FOUND);
1010 aperts = cpl_apertures_new_from_image(self, labels);
1013 code = irplib_closeset_aperture(aperts, x_initial, y_initial, &ifluxapert);
1016 cpl_apertures_delete(aperts);
1017 cpl_image_delete(labels);
1018 return cpl_error_set(cpl_func, CPL_ERROR_DATA_NOT_FOUND);
1021 npixobj = cpl_apertures_get_npix(aperts, ifluxapert);
1022 objradius = sqrt((
double)npixobj * CPL_MATH_1_PI);
1024 winsize = IRPLIB_MIN(IRPLIB_MIN(nx, ny), irplib_roundup_power2
1025 ((uint32_t)(3.0 * objradius + 0.5)));
1027 xposmax = cpl_apertures_get_maxpos_x(aperts, ifluxapert);
1028 yposmax = cpl_apertures_get_maxpos_y(aperts, ifluxapert);
1029 xposcen = cpl_apertures_get_centroid_x(aperts, ifluxapert);
1030 yposcen = cpl_apertures_get_centroid_y(aperts, ifluxapert);
1031 valmax = cpl_apertures_get_max(aperts, ifluxapert);
1033 cpl_apertures_delete(aperts);
1034 cpl_image_delete(labels);
1036 cpl_msg_debug(cpl_func,
"Object radius at S/R=%g: %g (window-size=%u)",
1037 sigma, objradius, (
unsigned)winsize);
1038 cpl_msg_debug(cpl_func,
"Object-peak @ (%d, %d) = %g", (
int)xposmax,
1039 (
int)yposmax, valmax);
1041 gauss_parameters = cpl_array_new(7, CPL_TYPE_DOUBLE);
1042 cpl_array_set_double(gauss_parameters, 0, median);
1044 code = cpl_fit_image_gaussian(self, NULL, xposcen, yposcen,
1045 winsize, winsize, gauss_parameters,
1050 const double M_x = cpl_array_get_double(gauss_parameters, 3, NULL);
1051 const double M_y = cpl_array_get_double(gauss_parameters, 4, NULL);
1053 valfit = irplib_gaussian_eval_2d(gauss_parameters, M_x, M_y);
1055 if (!cpl_errorstate_is_equal(prestate)) {
1056 code = cpl_error_get_code();
1062 cpl_msg_debug(cpl_func,
"Gauss-fit @ (%g, %g) = %g",
1066 cpl_array_delete(gauss_parameters);
1068#ifdef IRPLIB_STREHL_USE_CPL_IMAGE_FIT_GAUSSIAN
1069 if (code || valfit < valmax) {
1070 cpl_errorstate_set(prestate);
1072 code = cpl_image_fit_gaussian(self, xposcen, yposcen,
1073 (
int)(2.0 * objradius),
1083 valfit = irplib_gaussian_2d(0.0, 0.0, norm, sig_x, sig_y);
1085 cpl_msg_debug(cpl_func,
"Gauss-Fit @ (%g, %g) = %g. norm=%g, "
1086 "sigma=(%g, %g)", xcen, ycen, valfit, norm,
1089 if (valfit > valmax) {
1098 if (code || valfit < valmax) {
1099 cpl_errorstate_set(prestate);
1105 return code ? cpl_error_set_where(cpl_func) : CPL_ERROR_NONE;