28#if !defined(_XOPEN_SOURCE) || (_XOPEN_SOURCE - 0) < 600
29#define _XOPEN_SOURCE 600
32#include "hdrl_strehl.h"
33#include "hdrl_image.h"
34#include "hdrl_types.h"
35#include "hdrl_utils.h"
47static hdrl_strehl_result bad_result = {
95 double bkg_radius_low;
96 double bkg_radius_high;
97} hdrl_strehl_parameter;
100static hdrl_parameter_typeobj hdrl_strehl_parameter_type = {
101 HDRL_PARAMETER_STREHL,
102 (hdrl_alloc *)&cpl_malloc,
103 (hdrl_free *)&cpl_free,
105 sizeof(hdrl_strehl_parameter),
119hdrl_strehl_parameter_verify(
const hdrl_parameter * param)
121 cpl_error_ensure(param != NULL, CPL_ERROR_NULL_INPUT,
122 return CPL_ERROR_NULL_INPUT,
"NULL Input Parameters");
124 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
125 "Expected Strehl parameter") ;
127 const hdrl_strehl_parameter * param_loc = (
const hdrl_strehl_parameter *)param ;
129 cpl_error_ensure(param_loc->wavelength >= 0, CPL_ERROR_ILLEGAL_INPUT,
130 return CPL_ERROR_ILLEGAL_INPUT,
"wavelength must be >=0");
131 cpl_error_ensure(param_loc->m1 >= 0, CPL_ERROR_ILLEGAL_INPUT,
132 return CPL_ERROR_ILLEGAL_INPUT,
"m1 radius must be >=0");
133 cpl_error_ensure(param_loc->m2 >= 0, CPL_ERROR_ILLEGAL_INPUT,
134 return CPL_ERROR_ILLEGAL_INPUT,
"m2 radius must be >=0");
135 cpl_error_ensure(param_loc->m1 > param_loc->m2, CPL_ERROR_ILLEGAL_INPUT,
136 return CPL_ERROR_ILLEGAL_INPUT,
137 "m1 radius must be larger than m2 radius");
138 cpl_error_ensure(param_loc->pixel_scale_x >= 0, CPL_ERROR_ILLEGAL_INPUT,
139 return CPL_ERROR_ILLEGAL_INPUT,
"pixel_scale_x must be >=0");
140 cpl_error_ensure(param_loc->pixel_scale_y >= 0, CPL_ERROR_ILLEGAL_INPUT,
141 return CPL_ERROR_ILLEGAL_INPUT,
"pixel_scale_y must be >=0");
143 cpl_error_ensure(param_loc->flux_radius >= 0, CPL_ERROR_ILLEGAL_INPUT,
144 return CPL_ERROR_ILLEGAL_INPUT,
"flux_radius must be >=0");
146 cpl_error_ensure(param_loc->m1 >= param_loc->m2, CPL_ERROR_ILLEGAL_INPUT,
147 return CPL_ERROR_ILLEGAL_INPUT,
"m1 must be >=m2");
148 if (param_loc->bkg_radius_low > 0) {
149 cpl_error_ensure(param_loc->bkg_radius_low >= param_loc->flux_radius,
150 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
151 "bkg_radius_low must be >=flux_radius");
152 cpl_error_ensure(param_loc->bkg_radius_high > param_loc->bkg_radius_low,
153 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
154 "bkg_radius_high must be >bkg_radius_low");
157 cpl_error_ensure(param_loc->bkg_radius_high < 0,
158 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
159 "bkg_radius_high must be < 0 if bkg_radius_low is < 0");
162 return CPL_ERROR_NONE ;
189 double m1_radius,
double m2_radius,
190 double pixel_scale_x,
double pixel_scale_y,
double flux_radius,
191 double bkg_radius_low,
double bkg_radius_high) {
194 hdrl_strehl_parameter * p = (hdrl_strehl_parameter *)
195 hdrl_parameter_new(&hdrl_strehl_parameter_type);
197 p->wavelength = wavelength ;
200 p->pixel_scale_x = pixel_scale_x;
201 p->pixel_scale_y = pixel_scale_y;
202 p->flux_radius = flux_radius;
203 p->bkg_radius_low = bkg_radius_low;
204 p->bkg_radius_high = bkg_radius_high;
206 if (hdrl_strehl_parameter_verify((hdrl_parameter *)p)) {
210 return (hdrl_parameter *)p;
223 return hdrl_parameter_check_type(self, &hdrl_strehl_parameter_type);
234 const hdrl_parameter * p)
236 cpl_ensure(p, CPL_ERROR_NULL_INPUT, -1.0);
237 return p != NULL ? ((
const hdrl_strehl_parameter *)p)->wavelength : 0.;
248 const hdrl_parameter * p)
250 cpl_ensure(p, CPL_ERROR_NULL_INPUT, -1.0);
251 return p != NULL ? ((
const hdrl_strehl_parameter *)p)->m1 : 0.;
262 const hdrl_parameter * p)
264 cpl_ensure(p, CPL_ERROR_NULL_INPUT, -1.0);
265 return p != NULL ? ((
const hdrl_strehl_parameter *)p)->m2 : 0.;
276 const hdrl_parameter * p)
278 cpl_ensure(p, CPL_ERROR_NULL_INPUT, -1.0);
279 return p != NULL ? ((
const hdrl_strehl_parameter *)p)->pixel_scale_x : 0.;
290 const hdrl_parameter * p)
292 cpl_ensure(p, CPL_ERROR_NULL_INPUT, -1.0);
293 return p != NULL ? ((
const hdrl_strehl_parameter *)p)->pixel_scale_y : 0.;
304 const hdrl_parameter * p)
306 cpl_ensure(p, CPL_ERROR_NULL_INPUT, -1.0);
307 return p != NULL ? ((
const hdrl_strehl_parameter *)p)->flux_radius : 0.;
318 const hdrl_parameter * p)
320 cpl_ensure(p, CPL_ERROR_NULL_INPUT, -1.0);
321 return p != NULL ? ((
const hdrl_strehl_parameter *)p)->bkg_radius_low : 0.;
332 const hdrl_parameter * p)
334 cpl_ensure(p, CPL_ERROR_NULL_INPUT, -1.0);
335 return p != NULL ? ((
const hdrl_strehl_parameter *)p)->bkg_radius_high : 0.;
361 const char *base_context,
366 cpl_ensure(prefix && base_context && par,
367 CPL_ERROR_NULL_INPUT, NULL);
370 CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
372 cpl_parameterlist *parlist = cpl_parameterlist_new();
375 hdrl_setup_vparameter(parlist, prefix,
".",
"",
"wavelength", base_context,
376 "Wavelength [m].", CPL_TYPE_DOUBLE,
380 hdrl_setup_vparameter(parlist, prefix,
".",
"",
"m1", base_context,
381 "Telescope radius [m].", CPL_TYPE_DOUBLE,
385 hdrl_setup_vparameter(parlist, prefix,
".",
"",
"m2", base_context,
386 "Telescope obstruction radius [m].", CPL_TYPE_DOUBLE,
390 hdrl_setup_vparameter(parlist, prefix,
".",
"",
"pixel-scale-x", base_context,
391 "Detector X pixel scale on sky [arcsec].",
396 hdrl_setup_vparameter(parlist, prefix,
".",
"",
"pixel-scale-y", base_context,
397 "Detector Y pixel scale on sky [arcsec].",
402 hdrl_setup_vparameter(parlist, prefix,
".",
"",
"flux-radius", base_context,
403 "PSF Flux integration radius [arcsec].",
408 hdrl_setup_vparameter(parlist, prefix,
".",
"",
"bkg-radius-low", base_context,
409 "PSF background inner radii [arcsec].",
414 hdrl_setup_vparameter(parlist, prefix,
".",
"",
"bkg-radius-high", base_context,
415 "PSF background outer radius [arcsec].",
419 if (cpl_error_get_code()) {
420 cpl_parameterlist_delete(parlist);
449 const cpl_parameterlist * parlist,
452 cpl_ensure(prefix && parlist, CPL_ERROR_NULL_INPUT, NULL);
454 const cpl_parameter * par;
455 double wavelength, m1, m2, pixel_scale_x, pixel_scale_y,
456 flux_radius, bkg_radius_low,bkg_radius_high;
461 par=cpl_parameterlist_find_const(parlist, name);
462 wavelength = cpl_parameter_get_double(par);
467 par=cpl_parameterlist_find_const(parlist, name);
468 m1 = cpl_parameter_get_double(par);
473 par=cpl_parameterlist_find_const(parlist, name);
474 m2 = cpl_parameter_get_double(par);
479 par=cpl_parameterlist_find_const(parlist, name);
480 pixel_scale_x = cpl_parameter_get_double(par);
485 par=cpl_parameterlist_find_const(parlist, name);
486 pixel_scale_y = cpl_parameter_get_double(par);
491 par=cpl_parameterlist_find_const(parlist, name);
492 flux_radius = cpl_parameter_get_double(par);
497 par=cpl_parameterlist_find_const(parlist, name);
498 bkg_radius_low = cpl_parameter_get_double(par);
503 par=cpl_parameterlist_find_const(parlist, name);
504 bkg_radius_high = cpl_parameter_get_double(par);
508 if (cpl_error_get_code()) {
509 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
510 "Error while parsing parameterlist with prefix %s", prefix);
514 pixel_scale_x, pixel_scale_y,
515 flux_radius, bkg_radius_low, bkg_radius_high) ;
533hdrl_image_max_where(
const hdrl_image * himg, cpl_mask * mask)
555hdrl_image_sum_where(
const hdrl_image * himg, cpl_mask * mask)
575hdrl_image_median_where(
const hdrl_image * himg, cpl_mask * mask)
595hdrl_image_stdev_where(
const hdrl_image * himg, cpl_mask * mask)
602 return mad * CPL_MATH_STD_MAD;
640compute_psf(
double lam,
double m1,
double m2,
641 double pixscale_x,
double pixscale_y,
642 double cx,
double cy,
643 size_t nx,
size_t ny)
645 cpl_image * psf = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
646 double * data = cpl_image_get_data(psf);
648 double as_2_rad = CPL_MATH_2PI / (360. * 3600);
653 double centerx = (-(nx / 2.) + cx - 1 + 0.5) * pixscale_x;
654 double centery = (-(ny / 2.) + cy - 1 + 0.5) * pixscale_y;
655 double xhigh = ((nx - 1) * pixscale_x / 2) - centerx;
656 double yhigh = ((ny - 1) * pixscale_y / 2) - centery;
657 double xlow = -((nx - 1) * pixscale_x / 2) - centerx;
658 double ylow = -((ny - 1) * pixscale_y / 2) - centery;
659 double step_x = (xhigh - xlow) / (nx - 1);
660 double step_y = (yhigh - ylow) / (ny - 1);
662 HDRL_OMP(omp parallel
for)
663 for (
size_t iy = 0; iy < ny; iy++) {
664 double y = iy == ny - 1 ? yhigh : ylow + iy * step_y;
665 for (
size_t ix = 0; ix < nx; ix++) {
666 double x = ix == nx - 1 ? xhigh : xlow + ix * step_x;
667 double r = sqrt(x*x + y*y) * as_2_rad * CPL_MATH_2PI * m1 / lam;
669 data[iy * nx + ix] = 1.;
672 double airy = (2 * j1(r) / r - 2 * e * j1(e * r) / r);
673 double c = (1 - e * e);
674 data[iy * nx + ix] = 1 / (c * c) * airy * airy;
702static cpl_error_code apertures_find_max_flux(
const cpl_apertures * self,
703 int * ind,
int nfind)
705 const int nsize = cpl_apertures_get_size(self);
709 cpl_ensure_code(nsize > 0, cpl_error_get_code());
710 cpl_ensure_code(ind, CPL_ERROR_NULL_INPUT);
711 cpl_ensure_code(nfind > 0, CPL_ERROR_ILLEGAL_INPUT);
712 cpl_ensure_code(nfind <= nsize, CPL_ERROR_ILLEGAL_INPUT);
714 for (ifind=0; ifind < nfind; ifind++) {
718 for (i=1; i <= nsize; i++) {
722 for (k=0; k < ifind; k++)
if (ind[k] == i)
break;
726 const double flux = cpl_apertures_get_flux(self, i);
728 if (maxind < 0 || flux > maxflux) {
737 return CPL_ERROR_NONE;
761gaussian_maxpos(
const cpl_image * self,
768 const cpl_size nx = cpl_image_get_size_x(self);
769 const cpl_size ny = cpl_image_get_size_y(self);
773 const double median = cpl_image_get_median_dev(self, &med_dist);
774 cpl_mask * selection;
775 cpl_size nlabels = 0;
776 cpl_image * labels = NULL;
777 cpl_apertures * aperts;
781 cpl_size xposmax, yposmax;
782 double xposcen, yposcen;
783 double valmax, valfit = -1.0;
784 cpl_array * gauss_parameters = NULL;
785 cpl_errorstate prestate = cpl_errorstate_get();
786 cpl_error_code code = CPL_ERROR_NONE;
789 cpl_ensure_code( sigma > 0.0, CPL_ERROR_ILLEGAL_INPUT);
791 selection = cpl_mask_new(nx, ny);
794 for (; iretry > 0 && nlabels == 0; iretry--, sigma *= 0.5) {
797 const double threshold = median + sigma * med_dist;
801 code = cpl_mask_threshold_image(selection, self, threshold, DBL_MAX,
807 cpl_image_delete(labels);
808 labels = cpl_image_labelise_mask_create(selection, &nlabels);
812 cpl_mask_delete(selection);
815 cpl_image_delete(labels);
816 return cpl_error_set_where(cpl_func);
817 }
else if (nlabels == 0) {
818 cpl_image_delete(labels);
819 return cpl_error_set(cpl_func, CPL_ERROR_DATA_NOT_FOUND);
822 aperts = cpl_apertures_new_from_image(self, labels);
825 code = apertures_find_max_flux(aperts, &ifluxapert, 1);
828 cpl_apertures_delete(aperts);
829 cpl_image_delete(labels);
830 return cpl_error_set(cpl_func, CPL_ERROR_DATA_NOT_FOUND);
833 npixobj = cpl_apertures_get_npix(aperts, ifluxapert);
834 objradius = sqrt((
double)npixobj * CPL_MATH_1_PI);
835 winsize = CX_MIN(CX_MIN(nx, ny), (3.0 * objradius));
837 xposmax = cpl_apertures_get_maxpos_x(aperts, ifluxapert);
838 yposmax = cpl_apertures_get_maxpos_y(aperts, ifluxapert);
839 xposcen = cpl_apertures_get_centroid_x(aperts, ifluxapert);
840 yposcen = cpl_apertures_get_centroid_y(aperts, ifluxapert);
841 valmax = cpl_apertures_get_max(aperts, ifluxapert);
843 cpl_apertures_delete(aperts);
844 cpl_image_delete(labels);
846 cpl_msg_debug(cpl_func,
"Object radius at S/R=%g: %g (window-size=%u)",
847 sigma, objradius, (
unsigned)winsize);
848 cpl_msg_debug(cpl_func,
"Object-peak @ (%d, %d) = %g", (
int)xposmax,
849 (
int)yposmax, valmax);
853 gauss_parameters = cpl_array_new(7, CPL_TYPE_DOUBLE);
854 cpl_array_set_double(gauss_parameters, 0, median);
856 code = cpl_fit_image_gaussian(self, NULL, xposmax, yposmax,
857 winsize, winsize, gauss_parameters,
862 const double M_x = cpl_array_get_double(gauss_parameters, 3, NULL);
863 const double M_y = cpl_array_get_double(gauss_parameters, 4, NULL);
865 valfit = hcpl_gaussian_eval_2d(gauss_parameters, M_x, M_y);
867 if (!cpl_errorstate_is_equal(prestate)) {
868 code = cpl_error_get_code();
874 cpl_msg_debug(cpl_func,
"Gauss-fit @ (%g, %g) = %g",
878 cpl_array_delete(gauss_parameters);
880 if (code || valfit < valmax) {
881 cpl_errorstate_set(prestate);
887 return code ? cpl_error_set_where(cpl_func) : CPL_ERROR_NONE;
904strehl_disk_mask(
const cpl_image * im,
909 const intptr_t nx = cpl_image_get_size_x(im);
910 const intptr_t ny = cpl_image_get_size_y(im);
912 const intptr_t lx = (intptr_t)(xpos - rad);
913 const intptr_t ly = (intptr_t)(ypos - rad);
915 const intptr_t ux = (intptr_t)(xpos + rad) + 1;
916 const intptr_t uy = (intptr_t)(ypos + rad) + 1;
918 const double sqr = rad * rad;
923 cpl_ensure(im != NULL, CPL_ERROR_NULL_INPUT, NULL);
924 cpl_ensure(rad > 0.0, CPL_ERROR_ILLEGAL_INPUT, NULL);
926 m = cpl_mask_new(nx, ny);
928 for (intptr_t j = CX_MAX(ly, 0); j < CX_MIN(uy, ny); j++) {
929 const double yj = (double)j - ypos;
930 for (intptr_t i = CX_MAX(lx, 0); i < CX_MIN(ux, nx); i++) {
931 const double xi = (double)i - xpos;
932 const double dist = yj * yj + xi * xi;
934 if (!cpl_image_is_rejected(im, i + 1, j + 1)) {
935 cpl_mask_set(m, i + 1, j + 1, CPL_BINARY_1);
953static inline cpl_image * hdrl_rebin(cpl_image * img,
size_t sampling)
955 cpl_size lnx = cpl_image_get_size_x(img);
956 cpl_size lny = cpl_image_get_size_y(img);
957 cpl_size nx = lnx / sampling;
958 cpl_size ny = lny / sampling;
959 cpl_image * n = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
960 double * ld = cpl_image_get_data_double(img);
961 double * nd = cpl_image_get_data_double(n);
962 for (
size_t iy = 0; iy < (size_t)ny; iy++) {
963 for (
size_t ix = 0; ix < (size_t)nx; ix++) {
964 for (
size_t ly = 0; ly < sampling; ly++) {
965 for (
size_t lx = 0; lx < sampling; lx++) {
966 nd[iy * nx + ix] += ld[((iy*sampling) + ly) * lnx + lx + (ix*sampling)];
1002compute_psf_flux(
double lam,
double m1,
double m2,
double r)
1006 r *= CPL_MATH_2PI * m1 / lam;
1007 double E = 1 - j0(r) * j0(r) - j1(r) * j1(r) + e*e *
1008 (1 - j0(e*r)*j0(e*r) - j1(e*r)*j1(e*r));
1009 E *= 1 / (1 - e * e);
1074static hdrl_strehl_result
1075compute_strehl2(
const hdrl_image * himg,
double lam,
1076 double m1_radius,
double m2_radius,
1077 double pixscale_x,
double pixscale_y,
1078 double peak_x,
double peak_y,
1079 double radius_arcsec)
1082 double min_pscale = CX_MIN(pixscale_x, pixscale_y);
1083 double radius_pix = (radius_arcsec / min_pscale);
1086 cpl_size wins = 2 * radius_pix;
1087 cpl_msg_debug(cpl_func,
"strehl psf window size %d", (
int)wins);
1088 double smallx = peak_x - (floor(peak_x) - wins/2);
1089 double smally = peak_y - (floor(peak_y) - wins/2);
1091 size_t sampling = 16;
1092 intptr_t nnx = wins * sampling;
1093 intptr_t nny = wins * sampling;
1094 cpl_image * lpsf = compute_psf(lam, m1_radius, m2_radius,
1095 pixscale_x / sampling, pixscale_y / sampling,
1096 smallx*sampling,smally*sampling,
1101 cpl_image * epsf = cpl_image_extract(lpsf, 1 + o, 1 + o, nnx - o, nny - o);
1102 cpl_image * psf = hdrl_rebin(epsf, sampling);
1103 cpl_image_delete(epsf);
1104 cpl_image_delete(lpsf);
1106 cpl_image_divide_scalar(psf, cpl_image_get_max(psf)/cpl_image_get_max(img));
1107 cpl_msg_debug(cpl_func,
"position/peak of data: %g %g", peak_x, peak_y);
1110 double xposfit, yposfit, peak;
1111 gaussian_maxpos(psf, 5, &xposfit, &yposfit, &peak);
1112 cpl_msg_debug(cpl_func,
"position/peak of psf: %g %g", xposfit, yposfit);
1116 cpl_mask * im = strehl_disk_mask(img, peak_x, peak_y, radius_pix);
1117 hdrl_value ipeak = hdrl_image_max_where(himg, im);
1118 cpl_msg_debug(cpl_func,
"Computing flux on %d pixel radius, total pixels %ld",
1120 (
long)(cpl_mask_get_size_x(im) * cpl_mask_get_size_y(im) -
1121 cpl_mask_count(im)));
1122 hdrl_value iflux = hdrl_image_sum_where(himg, im);
1123 cpl_msg_debug(cpl_func,
"flux ring/total data: %g (%g) %g", iflux.data,
1124 iflux.error, cpl_image_get_flux(img));
1125 cpl_mask_delete(im);
1127 double ratio_img = ipeak.data / iflux.data;
1130 double ppeak = cpl_image_get_max(psf);
1132 strehl_disk_mask(psf, wins / 2 - 1, wins / 2 - 1, radius_pix);
1134 hdrl_value pflux = hdrl_image_sum_where(tmpimg, pm);
1136 cpl_msg_debug(cpl_func,
"flux ring/total psf: %g %g", pflux.data,
1137 cpl_image_get_flux(psf));
1138 cpl_mask_delete(pm);
1140 double ratio_psf = ppeak / pflux.data;
1142 cpl_msg_debug(cpl_func,
"data peak,flux,ratio: %g %g: %g",
1143 ipeak.data, iflux.data, ratio_img);
1144 cpl_msg_debug(cpl_func,
"psf peak,flux,ratio: %g %g: %g",
1145 ppeak, pflux.data, ratio_psf);
1147 double strehl = ratio_img / ratio_psf;
1148 double strehl_err = strehl *
1149 sqrt(ipeak.error * ipeak.error / (ipeak.data * ipeak.data) +
1150 iflux.error * iflux.error / (iflux.data * iflux.data));
1152 cpl_msg_debug(cpl_func,
"Strehl ratio %g +/- %g", strehl, strehl_err);
1153 cpl_image_delete(psf);
1155 hdrl_strehl_result r;
1156 r.strehl_value = (hdrl_value){strehl, strehl_err};
1157 r.star_peak = ipeak;
1158 r.star_flux = iflux;
1159 r.star_background = (hdrl_value){0., 0.};
1165 r.computed_background_error = -1.;
1166 r.nbackground_pixels = 0;
1209static hdrl_strehl_result
1210compute_strehl(
const hdrl_image * himg_,
double lam,
1211 double m1_radius,
double m2_radius,
1212 double pixscale_x,
double pixscale_y,
double flux_radius,
1213 double bkg_radius_low,
double bkg_radius_high)
1216 double xposfit, yposfit, peak;
1217 double min_pscale = CX_MIN(pixscale_x, pixscale_y);
1221 cpl_msg_warning(cpl_func,
1222 "%zu bad pixels in strehl input, interpolating.",
1229 if (gaussian_maxpos(img, 5, &xposfit, &yposfit, &peak) != CPL_ERROR_NONE) {
1233 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
1234 "detected peak of star smaller than zero, "
1235 "gaussian fit likely failed to fit the star");
1238 hdrl_value bkg = {0, 0};
1239 double comp_bkg_error = -1;
1240 size_t nbkg_pix = 0;
1241 if ((bkg_radius_low < 0 && bkg_radius_high >= 0) ||
1242 (bkg_radius_low >= 0 && bkg_radius_high < 0)) {
1243 cpl_error_set_message(cpl_func, CPL_ERROR_INCOMPATIBLE_INPUT,
1244 "background radius parameters must be larger "
1245 "zero or both negative");
1248 else if (bkg_radius_low >= 0 && bkg_radius_high >= 0) {
1249 if (bkg_radius_low >= bkg_radius_high) {
1250 cpl_error_set_message(cpl_func, CPL_ERROR_INCOMPATIBLE_INPUT,
1251 "low background radius parameters must be smaller "
1252 "than large background radius");
1256 cpl_mask * high = strehl_disk_mask(img, xposfit, yposfit,
1257 bkg_radius_high / min_pscale);
1258 cpl_mask * low = strehl_disk_mask(img, xposfit, yposfit,
1259 bkg_radius_low / min_pscale);
1260 cpl_mask_xor(low, high);
1261 nbkg_pix = cpl_mask_count(low);
1262 if (nbkg_pix == 0) {
1263 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
1264 "No valid pixels in background");
1265 cpl_mask_delete(low);
1266 cpl_mask_delete(high);
1271 bkg = hdrl_image_median_where(himg, low);
1272 comp_bkg_error = hdrl_image_stdev_where(himg, low) / sqrt(nbkg_pix);
1274 cpl_msg_debug(cpl_func,
"Median estimated background: %g +- %g "
1275 "(computed error %g)", bkg.data, bkg.error,
1277 cpl_mask_delete(low);
1278 cpl_mask_delete(high);
1283 hdrl_strehl_result r = compute_strehl2(himg, lam, m1_radius, m2_radius,
1284 pixscale_x, pixscale_y, xposfit, yposfit,
1287 r.star_background = bkg;
1290 r.computed_background_error = comp_bkg_error;
1291 r.nbackground_pixels = nbkg_pix;
1328 hdrl_strehl_result r;
1331 cpl_error_ensure(himg && params, CPL_ERROR_NULL_INPUT,
1334 if (hdrl_strehl_parameter_verify(params) != CPL_ERROR_NONE) {
1339 const hdrl_strehl_parameter * p_loc = (hdrl_strehl_parameter *)params ;
1341 r = compute_strehl(himg,p_loc->wavelength, p_loc->m1, p_loc->m2,
1342 p_loc->pixel_scale_x, p_loc->pixel_scale_y,
1343 p_loc->flux_radius, p_loc->bkg_radius_low,
1344 p_loc->bkg_radius_high);
hdrl_value hdrl_image_get_pixel(const hdrl_image *self, cpl_size xpos, cpl_size ypos, int *pis_rejected)
get pixel values of hdrl_image
cpl_error_code hdrl_image_reject_from_mask(hdrl_image *self, const cpl_mask *map)
set bpm of hdrl_image
hdrl_value hdrl_image_get_median(const hdrl_image *self)
computes the median and associated error of an image.
hdrl_image * hdrl_image_duplicate(const hdrl_image *himg)
copy hdrl_image
hdrl_value hdrl_image_get_sum(const hdrl_image *self)
computes the sum of all pixel values and the associated error of an image.
cpl_image * hdrl_image_get_error(hdrl_image *himg)
get error as cpl image
cpl_size hdrl_image_count_rejected(const hdrl_image *self)
return number of rejected pixels
hdrl_image * hdrl_image_create(const cpl_image *image, const cpl_image *error)
create a new hdrl_image from to existing images by copying them
cpl_image * hdrl_image_get_image(hdrl_image *himg)
get data as cpl image
const cpl_image * hdrl_image_get_image_const(const hdrl_image *himg)
get data as cpl image
void hdrl_image_delete(hdrl_image *himg)
delete hdrl_image
cpl_error_code hdrl_image_sub_scalar(hdrl_image *self, hdrl_value value)
Elementwise subtraction of a scalar from an image.
double hdrl_strehl_parameter_get_m1(const hdrl_parameter *p)
Access the primary mirror radius in the Strehl parameter.
cpl_parameterlist * hdrl_strehl_parameter_create_parlist(const char *base_context, const char *prefix, hdrl_parameter *par)
Create parameter list for the Strehl computation.
hdrl_strehl_result hdrl_strehl_compute(const hdrl_image *himg, hdrl_parameter *params)
This function computes the Strehl ratio.
double hdrl_strehl_parameter_get_flux_radius(const hdrl_parameter *p)
Access the total flux radius in the Strehl parameter.
double hdrl_strehl_parameter_get_wavelength(const hdrl_parameter *p)
Access the wavelength in the Strehl parameter.
double hdrl_strehl_parameter_get_m2(const hdrl_parameter *p)
Access the obstruction radius in the Strehl parameter.
double hdrl_strehl_parameter_get_bkg_radius_high(const hdrl_parameter *p)
Access the background region external radius in the Strehl parameter.
hdrl_parameter * hdrl_strehl_parameter_create(double wavelength, double m1_radius, double m2_radius, double pixel_scale_x, double pixel_scale_y, double flux_radius, double bkg_radius_low, double bkg_radius_high)
Creates Strehl Parameters object.
double hdrl_strehl_parameter_get_pixel_scale_x(const hdrl_parameter *p)
Access the image X pixel scale in the Strehl parameter.
double hdrl_strehl_parameter_get_bkg_radius_low(const hdrl_parameter *p)
Access the background region internal radius in the Strehl parameter.
hdrl_parameter * hdrl_strehl_parameter_parse_parlist(const cpl_parameterlist *parlist, const char *prefix)
Parse parameter list to create input parameters for the Strehl.
cpl_boolean hdrl_strehl_parameter_check(const hdrl_parameter *self)
Check that the parameter is a Strehl parameter.
double hdrl_strehl_parameter_get_pixel_scale_y(const hdrl_parameter *p)
Access the image Y pixel scale in the Strehl parameter.
char * hdrl_join_string(const char *sep_, int n,...)
join strings together