31#include <gsl/gsl_vector.h>
32#include <gsl/gsl_blas.h>
33#include <gsl/gsl_multifit_nlin.h>
35#include "moo_params.h"
36#include "moo_single.h"
37#include "moo_psf_single.h"
38#include "moo_badpix.h"
39#include "moo_fibres_table.h"
40#include "moo_model_flat.h"
56moo_image_get_linear_y_interpolated(cpl_image *image,
double y,
int x,
int *rej)
60 cpl_ensure(image != NULL, CPL_ERROR_NULL_INPUT, NAN);
61 cpl_ensure(rej != NULL, CPL_ERROR_NULL_INPUT, NAN);
63 int nx = cpl_image_get_size_x(image);
64 int ny = cpl_image_get_size_x(image);
65 cpl_ensure(x >= 1 && x <= nx, CPL_ERROR_ILLEGAL_INPUT, NAN);
66 cpl_ensure(y >= 1 && y <= ny, CPL_ERROR_ILLEGAL_INPUT, NAN);
68 int y1 = (int)floor(y);
71 double flux1 = cpl_image_get(image, x, y1, rej);
73 double flux2 = cpl_image_get(image, x, y2, rej);
75 if (*rej == CPL_BINARY_0) {
76 flux = flux1 + (y - y1) * (flux2 - flux1);
82moo_image_get_pixel_y_interpolated(cpl_image *image,
90 cpl_ensure(image != NULL, CPL_ERROR_NULL_INPUT, NAN);
91 cpl_ensure(rej != NULL, CPL_ERROR_NULL_INPUT, NAN);
93 int nx = cpl_image_get_size_x(image);
94 int ny = cpl_image_get_size_x(image);
95 cpl_ensure(x >= 1 && x <= nx, CPL_ERROR_ILLEGAL_INPUT, NAN);
96 cpl_ensure(y >= 1 && y <= ny, CPL_ERROR_ILLEGAL_INPUT, NAN);
98 double ystart = y - step / 2.0;
99 double ystop = y + step / 2.0;
101 int pix1 = (int)floor(ystart + 0.5);
102 int pix2 = (int)floor(ystop + 0.5);
105 flux = cpl_image_get(image, x, pix1, rej);
108 double frac1 = (pix1 + 0.5 - ystart) / step;
109 double frac2 = (ystop - (pix2 - 0.5)) / step;
110 double flux1 = cpl_image_get(image, x, pix1, rej);
111 double flux2 = cpl_image_get(image, x, pix2, rej);
113 flux = flux1 * frac1 + flux2 * frac2;
117 if (*rej == CPL_BINARY_1) {
125moo_image_get_pixel2_y_interpolated(cpl_image *image,
double y,
int x,
int *rej)
129 cpl_ensure(image != NULL, CPL_ERROR_NULL_INPUT, NAN);
130 cpl_ensure(rej != NULL, CPL_ERROR_NULL_INPUT, NAN);
132 int nx = cpl_image_get_size_x(image);
133 int ny = cpl_image_get_size_x(image);
134 cpl_ensure(x >= 1 && x <= nx, CPL_ERROR_ILLEGAL_INPUT, NAN);
135 cpl_ensure(y >= 1 && y <= ny, CPL_ERROR_ILLEGAL_INPUT, NAN);
137 int pix = (int)floor(y + 0.5);
139 flux = cpl_image_get(image, x, pix, rej);
141 if (*rej == CPL_BINARY_1) {
149moo_image_get_pixel_yerror_interpolated(cpl_image *error,
double y,
int x)
153 cpl_ensure(error != NULL, CPL_ERROR_NULL_INPUT, NAN);
155 int nx = cpl_image_get_size_x(error);
156 int ny = cpl_image_get_size_x(error);
157 cpl_ensure(x >= 1 && x <= nx, CPL_ERROR_ILLEGAL_INPUT, NAN);
158 cpl_ensure(y >= 1 && y <= ny, CPL_ERROR_ILLEGAL_INPUT, NAN);
160 int y1 = (int)floor(y + 0.5);
163 err = cpl_image_get(error, x, y1, &rej);
169moo_image_get_linear_yerror_interpolated(cpl_image *error,
double y,
int x)
173 cpl_ensure(error != NULL, CPL_ERROR_NULL_INPUT, NAN);
175 int nx = cpl_image_get_size_x(error);
176 int ny = cpl_image_get_size_x(error);
177 cpl_ensure(x >= 1 && x <= nx, CPL_ERROR_ILLEGAL_INPUT, NAN);
178 cpl_ensure(y >= 1 && y <= ny, CPL_ERROR_ILLEGAL_INPUT, NAN);
180 int y1 = (int)floor(y);
184 double err1 = cpl_image_get(error, x, y1, &rej);
186 double err2 = cpl_image_get(error, x, y2, &rej);
188 err = sqrt((y - y1) * err1 * (y - y1) * err1 +
189 (y2 - y) * err2 * (y2 - y) * err2);
195_moo_fit_model_flat_fibre(cpl_image *rec_fluxes,
int deg_poly)
197 int nx = cpl_image_get_size_x(rec_fluxes);
198 int ny = cpl_image_get_size_y(rec_fluxes);
200 cpl_bivector *points = cpl_bivector_new(nx);
201 cpl_vector *xpos = cpl_bivector_get_x(points);
202 cpl_vector *fluxes = cpl_bivector_get_y(points);
204 int *flag = cpl_calloc(nx,
sizeof(
int));
206 for (
int x = 0; x < nx; x++) {
207 cpl_vector_set(xpos, x, x + 1);
210 double *data = cpl_image_get_data(rec_fluxes);
211 for (
int y = 0; y < ny; y++) {
212 for (
int x = 0; x < nx; x++) {
213 double val = data[x + y * nx];
214 cpl_vector_set(fluxes, x, val);
227 for (
int x = 0; x < nx; x++) {
228 double val = cpl_vector_get(fluxes, x);
229 data[x + y * nx] = val;
233 cpl_bivector_delete(points);
238_moo_smooth_model_flat_fibre(cpl_image *rec_fluxes,
double s)
240 int nx = cpl_image_get_size_x(rec_fluxes);
241 int ny = cpl_image_get_size_y(rec_fluxes);
243 cpl_vector *fluxes = cpl_vector_new(nx);
245 double *data = cpl_image_get_data(rec_fluxes);
246 for (
int y = 0; y < ny; y++) {
248 for (
int x = 0; x < nx; x++) {
249 double val = data[x + y * nx];
253 cpl_vector_set(fluxes, x, val);
256 cpl_vector *filter = cpl_vector_new(nx - nbnan);
259 for (
int x = 0; x < nx; x++) {
260 double val = cpl_vector_get(fluxes, x);
263 cpl_vector_set(filter, filteridx, val);
268 cpl_vector *res = moo_hpfilter(filter, s);
272 for (
int x = 0; x < nx; x++) {
273 double val = cpl_vector_get(fluxes, x);
275 val = cpl_vector_get(res, x - nbnan);
280 data[x + y * nx] = val;
282 cpl_vector_delete(filter);
283 cpl_vector_delete(res);
286 cpl_vector_delete(fluxes);
290_moo_median_filter_model_flat_fibre(cpl_image *rec_fluxes,
294 int nx = cpl_image_get_size_x(rec_fluxes);
295 int ny = cpl_image_get_size_y(rec_fluxes);
297 cpl_mask *mask = cpl_mask_new(nx, ny);
299 cpl_mask *orig = cpl_image_get_bpm(rec_fluxes);
301 cpl_vector *fluxes = cpl_vector_new(nx);
302 int window_size = winhsize * 2 + 1;
304 cpl_vector *winflux = cpl_vector_new(window_size);
305 for (
int y = 1; y <= ny; y++) {
306 for (
int x = 1; x <= nx; x++) {
308 double val = cpl_image_get(rec_fluxes, x, y, &rej);
309 cpl_vector_set(fluxes, x - 1, val);
312 for (
int x = 0; x < window_size; x++) {
313 double val = cpl_vector_get(fluxes, x);
314 cpl_vector_set(winflux, x, val);
319 double sigma = (p2 - p1) / 2;
321 for (
int x = 0; x < winhsize + 1; x++) {
322 double val = cpl_vector_get(fluxes, x);
323 double diff = val - median;
324 if (fabs(diff) > kappa * sigma) {
325 cpl_mask_set(mask, x + 1, y, CPL_BINARY_1);
328 for (
int x = winhsize + 2; x < nx - winhsize - 1; x++) {
329 for (
int i = x - winhsize; i <= x + winhsize; i++) {
330 double val = cpl_vector_get(fluxes, i);
331 cpl_vector_set(winflux, i - x + winhsize, val);
333 double val = cpl_vector_get(winflux, winhsize);
338 sigma = (p2 - p1) / 2;
340 double diff = val - median;
341 if (fabs(diff) > kappa * sigma) {
342 cpl_mask_set(mask, x + 1, y, CPL_BINARY_1);
347 int count = cpl_mask_count(mask);
348 cpl_msg_info(__func__,
"new cosmics %d", count);
349 cpl_mask_or(orig, mask);
350 cpl_mask_delete(mask);
351 cpl_vector_delete(fluxes);
352 cpl_vector_delete(winflux);
356_moo_median_filter2_model_flat_fibre(cpl_image *rec_fluxes,
int winhsize)
358 int nx = cpl_image_get_size_x(rec_fluxes);
359 int ny = cpl_image_get_size_y(rec_fluxes);
361 cpl_vector *fluxes = cpl_vector_new(nx);
363 for (
int y = 1; y <= ny; y++) {
364 for (
int x = 1; x <= nx; x++) {
366 double val = cpl_image_get(rec_fluxes, x, y, &rej);
367 cpl_vector_set(fluxes, x - 1, val);
372 for (
int x = 1; x <= nx; x++) {
373 double val = cpl_vector_get(res, x - 1);
374 cpl_image_set(rec_fluxes, x, y, val);
377 cpl_vector_delete(res);
379 cpl_vector_delete(fluxes);
383_moo_savgol_model_flat_fibre(cpl_image *rec_fluxes,
387 int nx = cpl_image_get_size_x(rec_fluxes);
388 int ny = cpl_image_get_size_y(rec_fluxes);
390 cpl_vector *fluxes = cpl_vector_new(nx);
392 for (
int y = 1; y <= ny; y++) {
393 for (
int x = 1; x <= nx; x++) {
395 double val = cpl_image_get(rec_fluxes, x, y, &rej);
400 cpl_vector_set(fluxes, x - 1, val);
403 char *testname = cpl_sprintf(
"befsavgol.txt");
404 FILE *test = fopen(testname,
"w");
405 fprintf(test,
"#x f\n");
407 for (
int x = 0; x < nx; x++) {
408 double val = cpl_vector_get(fluxes, x);
409 fprintf(test,
"%d %f\n", x, val);
417 char *testname = cpl_sprintf(
"savgol.txt");
418 FILE *test = fopen(testname,
"w");
419 fprintf(test,
"#x f\n");
421 for (
int x = 0; x < nx; x++) {
422 double val = cpl_vector_get(res, x);
423 fprintf(test,
"%d %f\n", x, val);
430 for (
int x = 0; x < nx; x++) {
431 double val = cpl_vector_get(res, x);
432 cpl_image_set(rec_fluxes, x + 1, y, val);
434 cpl_vector_delete(res);
437 cpl_vector_delete(fluxes);
441_moo_extract_error_model_flat_fibre(hdrl_image *image,
447 cpl_image *rec_fluxes)
449 int nx = hdrl_image_get_size_x(image);
450 cpl_image *error = hdrl_image_get_error(image);
452 for (
int x = 1; x <= nx; x++) {
453 double centroid = centroids[(numfib - 1) * nx + x - 1];
455 if (!isnan(centroid)) {
456 for (
int y = -nstep_lo; y <= nstep_up; y++) {
457 double yc = centroid + y * step;
459 moo_image_get_pixel_yerror_interpolated(error, yc, x);
460 cpl_image_set(rec_fluxes, x, y + nstep_lo + 1, flux);
464 for (
int y = -nstep_lo; y <= nstep_up; y++) {
465 cpl_image_set(rec_fluxes, x, y + nstep_lo + 1, NAN);
472_moo_rectify_fibre(hdrl_image *image,
478 cpl_image *rec_fluxes)
480 int nx = hdrl_image_get_size_x(image);
481 cpl_image *data = hdrl_image_get_image(image);
483 for (
int x = 1; x <= nx; x++) {
484 double centroid = centroids[(numfib - 1) * nx + x - 1];
486 if (!isnan(centroid)) {
487 for (
int y = -nstep_lo; y <= nstep_up; y++) {
488 double yc = centroid + y * step;
493 moo_image_get_pixel2_y_interpolated(data, yc, x, &rej);
494 cpl_image_set(rec_fluxes, x, y + nstep_lo + 1, flux);
498 for (
int y = -nstep_lo; y <= nstep_up; y++) {
499 cpl_image_set(rec_fluxes, x, y + nstep_lo + 1, NAN);
506_moo_model_f(
const gsl_vector *px,
void *data, gsl_vector *f)
508 cpl_bivector *bvdata = (cpl_bivector *)data;
509 size_t n = cpl_bivector_get_size(bvdata);
510 cpl_vector *vy = cpl_bivector_get_y(bvdata);
511 cpl_vector *vx = cpl_bivector_get_x(bvdata);
513 double sigma = gsl_vector_get(px, 0);
514 double cexp = gsl_vector_get(px, 1);
515 double A = gsl_vector_get(px, 2);
516 double c0 = gsl_vector_get(px, 3);
517 double b = gsl_vector_get(px, 4);
519 moo_model *m = moo_model_new(sigma, cexp, A, c0, b);
523 for (i = 0; i < n; i++) {
525 double x = cpl_vector_get(vx, i);
526 double y = cpl_vector_get(vy, i);
528 double Yx = moo_model_eval(m, x);
530 gsl_vector_set(f, i, Yx - y);
538_moo_model_df(
const gsl_vector *px,
void *data, gsl_matrix *J)
540 cpl_bivector *bvdata = (cpl_bivector *)data;
541 size_t n = cpl_bivector_get_size(bvdata);
542 cpl_vector *vx = cpl_bivector_get_x(bvdata);
544 double sigma = gsl_vector_get(px, 0);
545 double cexp = gsl_vector_get(px, 1);
546 double A = gsl_vector_get(px, 2);
547 double c0 = gsl_vector_get(px, 3);
548 double b = gsl_vector_get(px, 4);
550 moo_model *m = moo_model_new(sigma, cexp, A, c0, b);
553 for (i = 0; i < n; i++) {
558 double x = cpl_vector_get(vx, i);
559 double dA = moo_model_dA_eval(m, x);
560 double dcexp = moo_model_dcexp_eval(m, x);
561 double dsigma = moo_model_dsigma_eval(m, x);
562 double dc0 = moo_model_dc0_eval(m, x);
564 gsl_matrix_set(J, i, 0, dsigma);
565 gsl_matrix_set(J, i, 1, dcexp);
566 gsl_matrix_set(J, i, 2, dA);
567 gsl_matrix_set(J, i, 3, dc0);
568 gsl_matrix_set(J, i, 4, 1.0);
576_moo_model_df2(
const gsl_vector *px,
void *data, gsl_matrix *J)
578 cpl_bivector *bvdata = (cpl_bivector *)data;
579 size_t n = cpl_bivector_get_size(bvdata);
580 cpl_vector *vx = cpl_bivector_get_x(bvdata);
582 double sigma = gsl_vector_get(px, 0);
583 double cexp = gsl_vector_get(px, 1);
584 double A = gsl_vector_get(px, 2);
585 double c0 = gsl_vector_get(px, 3);
586 double b = gsl_vector_get(px, 4);
588 moo_model *m = moo_model_new(sigma, cexp, A, c0, b);
591 for (i = 0; i < n; i++) {
596 double x = cpl_vector_get(vx, i);
597 double dA = moo_model_dA_eval(m, x);
598 double dc0 = moo_model_dc0_eval(m, x);
600 gsl_matrix_set(J, i, 0, 0);
601 gsl_matrix_set(J, i, 1, 0);
602 gsl_matrix_set(J, i, 2, dA);
603 gsl_matrix_set(J, i, 3, dc0);
604 gsl_matrix_set(J, i, 4, 1.0);
612_moo_model_fit(cpl_bivector *data,
613 int (*fit_func)(
const gsl_vector *,
void *, gsl_matrix *),
617 moo_model *model = NULL;
620 int N = cpl_bivector_get_size(data);
621 cpl_vector *dy = cpl_bivector_get_y(data);
622 double min = cpl_vector_get_min(dy);
623 double max = cpl_vector_get_max(dy);
624 double A = max - min;
626 const gsl_multifit_fdfsolver_type *T = gsl_multifit_fdfsolver_lmsder;
627 gsl_multifit_fdfsolver *s;
632 gsl_matrix *J = gsl_matrix_alloc(N, p);
633 gsl_matrix *covar = gsl_matrix_alloc(p, p);
635 gsl_multifit_function_fdf f;
636 double x_init[5] = { isigma, icexp, A, 0.0, b };
637 gsl_vector_view x = gsl_vector_view_array(x_init, p);
639 const double xtol = 1e-9;
640 const double gtol = 1e-9;
641 const double ftol = 0.0;
649 s = gsl_multifit_fdfsolver_alloc(T, N, p);
652 gsl_multifit_fdfsolver_set(s, &f, &x.vector);
655 gsl_multifit_fdfsolver_residual(s);
658 gsl_multifit_fdfsolver_driver(s, 50, xtol, gtol, ftol, &info);
660 gsl_multifit_fdfsolver_jac(s, J);
661 gsl_multifit_covar(J, 0.0, covar);
683 double sigma = gsl_vector_get(s->x, 0);
684 double cexp = gsl_vector_get(s->x, 1);
685 A = gsl_vector_get(s->x, 2);
686 double c0 = gsl_vector_get(s->x, 3);
687 b = gsl_vector_get(s->x, 4);
689 model = moo_model_new(sigma, cexp, A, c0, b);
690 gsl_multifit_fdfsolver_free(s);
691 gsl_matrix_free(covar);
698moo_model_new(
double sigma,
double cexp,
double A,
double c0,
double b)
700 moo_model *res = cpl_calloc(1,
sizeof(moo_model));
712moo_model_delete(moo_model *self)
720moo_model_eval(moo_model *model,
double x)
724 cpl_ensure(model != NULL, CPL_ERROR_NULL_INPUT, res);
726 double sigma = model->sigma;
727 double cexp = model->cexp;
729 double c0 = model->c0;
732 res = (A / (sqrt(CPL_MATH_2PI) * sigma)) *
733 exp(-0.5 * pow((fabs(x - c0) / sigma), cexp)) +
740moo_model_dA_eval(moo_model *model,
double x)
744 cpl_ensure(model != NULL, CPL_ERROR_NULL_INPUT, res);
746 double sigma = model->sigma;
747 double cexp = model->cexp;
748 double c0 = model->c0;
750 double absxc0 = fabs(x - c0);
752 res = (1 / (sqrt(CPL_MATH_2PI) * sigma)) *
753 exp(-0.5 * pow((absxc0 / sigma), cexp));
759moo_model_dsigma_eval(moo_model *model,
double x)
763 cpl_ensure(model != NULL, CPL_ERROR_NULL_INPUT, res);
765 double sigma = model->sigma;
766 double cexp = model->cexp;
767 double c0 = model->c0;
770 double absxc0 = fabs(x - c0);
772 res = -(A * pow(sigma, -cexp - 2) *
773 (2 * pow(sigma, cexp) - cexp * pow(absxc0, cexp)) *
774 exp(-pow(absxc0, cexp) / (2 * pow(sigma, cexp)))) /
775 (2 * sqrt(CPL_MATH_2PI));
780moo_model_dcexp_eval(moo_model *model,
double x)
784 cpl_ensure(model != NULL, CPL_ERROR_NULL_INPUT, res);
786 double sigma = model->sigma;
787 double cexp = model->cexp;
788 double c0 = model->c0;
795 double absxc0 = fabs(x - c0);
797 res = -(A * pow(sigma, -cexp - 1) * pow(absxc0, cexp) *
798 log(absxc0 / sigma) *
799 exp(-pow(absxc0, cexp) / (2. * pow(sigma, cexp)))) /
800 (2. * sqrt(CPL_MATH_2PI));
806moo_model_dc0_eval(moo_model *model,
double x)
810 cpl_ensure(model != NULL, CPL_ERROR_NULL_INPUT, res);
812 double sigma = model->sigma;
813 double cexp = model->cexp;
814 double c0 = model->c0;
817 double absxc0 = fabs(x - c0);
823 res = (A * cexp * pow(sigma, -cexp - 1) *
824 exp(-pow(absxc0, cexp) / (2 * pow(sigma, cexp))) *
826 (2 * sqrt(CPL_MATH_2PI) * (x - c0));
832_moo_model_fibre(hdrl_image *image,
840 cpl_image *rec_fluxes,
844 int nx = hdrl_image_get_size_x(image);
846 cpl_vector *vX = cpl_vector_new(nx);
847 cpl_vector *vCexp = cpl_vector_new(nx);
848 cpl_vector *vSigma = cpl_vector_new(nx);
852 for (
int x = win_hsize + 1; x <= nx - win_hsize; x += xstep) {
853 cpl_bivector *data = NULL;
855 for (
int x2 = x - win_hsize; x2 <= x + win_hsize; x2++) {
857 double centroid = centroids[(numfib - 1) * nx + x2 - 1];
858 double wlo = cpl_image_get(fwlo, x2, numfib, &prej);
859 double wup = cpl_image_get(fwup, x2, numfib, &prej);
860 int yl = floor(centroid - wlo);
861 int yu = ceil(centroid + wup);
863 int ysize = yu - yl + 1;
864 cpl_bivector *tmp_data = cpl_bivector_new(ysize);
865 cpl_vector *tmp_y = cpl_bivector_get_x(tmp_data);
866 cpl_vector *tmp_flux = cpl_bivector_get_y(tmp_data);
870 if (!isnan(centroid)) {
871 for (
int y = yl; y <= yu; y++) {
873 hdrl_value flux = hdrl_image_get_pixel(image, x2, y, &rej);
875 cpl_vector_set(tmp_y, tmp_size, y);
876 cpl_vector_set(tmp_flux, tmp_size, flux.data);
885 double frac = 1 - (double)nbrej / (
double)(tmp_size + nbrej);
887 if (tmp_size > 4 && frac > 0.8) {
888 cpl_vector_set_size(tmp_y, tmp_size);
889 cpl_vector_set_size(tmp_flux, tmp_size);
890 double min = cpl_vector_get_min(tmp_flux);
891 cpl_vector_subtract_scalar(tmp_y, centroid);
893 cpl_vector_subtract_scalar(tmp_flux, min);
895 double sum = cpl_vector_get_sum(tmp_flux);
898 cpl_vector_divide_scalar(tmp_flux, sum);
901 data = cpl_bivector_duplicate(tmp_data);
904 int size = cpl_bivector_get_size(data);
905 cpl_vector *dx = cpl_bivector_get_x(data);
906 cpl_vector *dy = cpl_bivector_get_y(data);
907 cpl_vector_set_size(dx, size + tmp_size);
908 cpl_vector_set_size(dy, size + tmp_size);
909 for (
int i = size; i < size + tmp_size; i++) {
910 cpl_vector_set(dx, i,
911 cpl_vector_get(tmp_y, i - size));
912 cpl_vector_set(dy, i,
913 cpl_vector_get(tmp_flux, i - size));
918 cpl_bivector_delete(tmp_data);
921 moo_model *model = _moo_model_fit(data, &_moo_model_df, 1., 2.);
924 cpl_vector_set(vX, model_size, x);
925 cpl_vector_set(vCexp, model_size, model->cexp);
926 cpl_vector_set(vSigma, model_size, model->sigma);
930 moo_model_delete(model);
931 cpl_bivector_delete(data);
934 cpl_vector_set_size(vX, model_size);
935 cpl_vector_set_size(vCexp, model_size);
936 cpl_vector_set_size(vSigma, model_size);
941 for (
int x = 1; x <= nx; x++) {
943 double centroid = centroids[(numfib - 1) * nx + x - 1];
944 double wlo = cpl_image_get(fwlo, x, numfib, &prej);
945 double wup = cpl_image_get(fwup, x, numfib, &prej);
946 int yl = floor(centroid - wlo);
947 int yu = ceil(centroid + wup);
949 int ysize = yu - yl + 1;
951 if (!isnan(centroid)) {
952 cpl_bivector *tmp_data = cpl_bivector_new(ysize);
953 cpl_vector *tmp_y = cpl_bivector_get_x(tmp_data);
954 cpl_vector *tmp_flux = cpl_bivector_get_y(tmp_data);
959 for (
int y = yl; y <= yu; y++) {
961 hdrl_value flux = hdrl_image_get_pixel(image, x, y, &rej);
963 cpl_vector_set(tmp_y, tmp_size, y);
964 cpl_vector_set(tmp_flux, tmp_size, flux.data);
971 double frac = (double)nbrej / (
double)(nbrej + tmp_size);
973 if (tmp_size > 4 && frac < 0.2) {
977 cpl_vector_set_size(tmp_flux, tmp_size);
978 cpl_vector_set_size(tmp_y, tmp_size);
979 cpl_vector_subtract_scalar(tmp_y, centroid);
981 int idx = cpl_vector_find(vX, x);
982 sigma = cpl_vector_get(resSigma, idx);
983 cexp = cpl_vector_get(resCexp, idx);
986 _moo_model_fit(tmp_data, &_moo_model_df2, sigma, cexp);
987 for (
int i = -nstep_lo; i <= nstep_up; i++) {
989 double v = moo_model_eval(model, y);
991 cpl_image_set(rec_fluxes, x, i + nstep_lo + 1, v);
994 moo_model_delete(model);
997 for (
int i = -nstep_lo; i <= nstep_up; i++) {
998 cpl_image_set(rec_fluxes, x, i + nstep_lo + 1, NAN);
1001 cpl_bivector_delete(tmp_data);
1004 for (
int i = -nstep_lo; i <= nstep_up; i++) {
1005 cpl_image_set(rec_fluxes, x, i + nstep_lo + 1, NAN);
1010 cpl_vector_delete(vX);
1011 cpl_vector_delete(vCexp);
1012 cpl_vector_delete(vSigma);
1014 cpl_vector_delete(resCexp);
1015 cpl_vector_delete(resSigma);
1018static moo_psf_single *
1019_moo_model_flat_single(moo_single *det,
1020 moo_loc_single *loc,
1021 const char *filename,
1022 moo_model_flat_params *params,
1026 cpl_errorstate prestate = cpl_errorstate_get();
1028 cpl_ensure(filename != NULL, CPL_ERROR_NULL_INPUT, NULL);
1029 cpl_ensure(det != NULL, CPL_ERROR_NULL_INPUT, NULL);
1036 int nb_fibres = cpl_image_get_size_y(fcentroids);
1037 int nx = cpl_image_get_size_x(fcentroids);
1039 cpl_image_reject_value(fwlo, CPL_VALUE_NAN);
1040 cpl_image_reject_value(fwup, CPL_VALUE_NAN);
1042 double max_fwlo = cpl_image_get_max(fwlo);
1043 double max_fwup = cpl_image_get_max(fwup);
1045 cpl_msg_indent_more();
1046 cpl_msg_info(
"moo_model_flat",
"Use Y localisation limits: %.3f,%.3f pix",
1047 -max_fwlo, max_fwup);
1049 int nstep_lo = (int)ceil(max_fwlo * params->oversamp);
1050 int nstep_up = (int)ceil(max_fwup * params->oversamp);
1052 int ny = nstep_lo + nstep_up + 1;
1056 double *centroids = cpl_image_get_data(fcentroids);
1057 double step = 1.0 / params->oversamp;
1059 cpl_image **imgtab = cpl_calloc(nb_fibres,
sizeof(cpl_image *));
1062 int win_hsize = params->winhsize;
1063 int xstep = params->xstep;
1065 cpl_msg_info(
"moo_model_flat",
1066 "Fitting profile shape parameters using step %d and window "
1067 "half size %d along X axis",
1069 cpl_msg_indent_less();
1071#if MOO_DEBUG_MODEL_FLAT
1072 int test_fibers[] = { 10, 250, 500 };
1073 for (
int ti = 306; ti <= 308; ti++) {
1077 cpl_sprintf(
"%s_%d_localise_step1.txt",
1080 cpl_sprintf(
"%s_%d_localise_align_centroid_step1.txt",
1082 char *testmodelname =
1083 cpl_sprintf(
"%s_%d_model_parameters_step1.txt",
1085 char *testmodel2name =
1086 cpl_sprintf(
"%s_%d_model_step1.txt",
1088 char *testmodel3name =
1089 cpl_sprintf(
"%s_%d_model_step2.txt",
1091 char *testmodelp2name =
1092 cpl_sprintf(
"%s_%d_model_parameters_step2.txt",
1094 FILE *test = fopen(testname,
"w");
1095 fprintf(test,
"#x x2 y flux err qual\n");
1097 FILE *test1 = fopen(testname1b,
"w");
1098 fprintf(test1,
"#x x2 y flux centroid min sum\n");
1100 FILE *test2 = fopen(testmodelname,
"w");
1101 fprintf(test2,
"#x centroid A sigma c0 cexp b\n");
1103 FILE *test3 = fopen(testmodel2name,
"w");
1104 fprintf(test3,
"#x y flux\n");
1106 FILE *test4 = fopen(testmodel3name,
"w");
1107 fprintf(test4,
"#x y yr flux\n");
1109 FILE *testp2 = fopen(testmodelp2name,
"w");
1110 fprintf(testp2,
"#x centroid A sigma c0 cexp b\n");
1112 cpl_vector *vX = cpl_vector_new(nx);
1113 cpl_vector *vA = cpl_vector_new(nx);
1114 cpl_vector *vCexp = cpl_vector_new(nx);
1115 cpl_vector *vSigma = cpl_vector_new(nx);
1116 cpl_vector *vB = cpl_vector_new(nx);
1117 cpl_vector *vC0 = cpl_vector_new(nx);
1121 for (
int x = win_hsize + 1; x <= nx - win_hsize; x += xstep) {
1122 cpl_bivector *data = NULL;
1124 for (
int x2 = x - win_hsize; x2 <= x + win_hsize; x2++) {
1126 double centroid = centroids[(numfib - 1) * nx + x2 - 1];
1127 double wlo = cpl_image_get(fwlo, x2, numfib, &prej);
1128 double wup = cpl_image_get(fwup, x2, numfib, &prej);
1129 int yl = floor(centroid - wlo);
1130 int yu = ceil(centroid + wup);
1132 int ysize = yu - yl + 1;
1133 cpl_bivector *tmp_data = cpl_bivector_new(ysize);
1134 cpl_vector *tmp_y = cpl_bivector_get_x(tmp_data);
1135 cpl_vector *tmp_flux = cpl_bivector_get_y(tmp_data);
1140 if (!isnan(centroid)) {
1141 for (
int y = yl; y <= yu; y++) {
1144 hdrl_image_get_pixel(image, x2, y, &rej);
1146 fprintf(test,
"%d %d %d %f %f %d\n", x, x2, y,
1147 flux.data, flux.error, rej);
1149 cpl_vector_set(tmp_y, tmp_size, y);
1150 cpl_vector_set(tmp_flux, tmp_size, flux.data);
1158 double frac = 1 - (double)nbrej / (
double)(tmp_size + nbrej);
1160 if (tmp_size >= 4 && frac > 0.9) {
1161 cpl_vector_set_size(tmp_y, tmp_size);
1162 cpl_vector_set_size(tmp_flux, tmp_size);
1164 double min = cpl_vector_get_min(tmp_flux);
1166 cpl_vector_subtract_scalar(tmp_y, centroid);
1167 cpl_vector_subtract_scalar(tmp_flux, min);
1168 double sum = cpl_vector_get_sum(tmp_flux);
1171 cpl_vector_divide_scalar(tmp_flux, sum);
1173 for (
int y = 0; y < tmp_size; y++) {
1174 double ry = cpl_vector_get(tmp_y, y);
1175 double rf = cpl_vector_get(tmp_flux, y);
1176 fprintf(test1,
"%d %d %f %f %f %f %f\n", x, x2, ry,
1177 rf, centroid, min, sum);
1181 data = cpl_bivector_duplicate(tmp_data);
1184 int size = cpl_bivector_get_size(data);
1185 cpl_vector *dx = cpl_bivector_get_x(data);
1186 cpl_vector *dy = cpl_bivector_get_y(data);
1187 cpl_vector_set_size(dx, size + tmp_size);
1188 cpl_vector_set_size(dy, size + tmp_size);
1189 for (
int i = size; i < size + tmp_size; i++) {
1190 cpl_vector_set(dx, i,
1191 cpl_vector_get(tmp_y, i - size));
1192 cpl_vector_set(dy, i,
1193 cpl_vector_get(tmp_flux,
1200 cpl_bivector_delete(tmp_data);
1204 moo_model *model = _moo_model_fit(data, &_moo_model_df, 1., 2.);
1205 if (model != NULL) {
1206 double centroid = centroids[(numfib - 1) * nx + x - 1];
1207 fprintf(test2,
"%d %f %f %f %f %f %f\n", x, centroid,
1208 model->A, model->sigma, model->c0, model->cexp,
1210 cpl_vector_set(vX, model_size, x);
1211 cpl_vector_set(vA, model_size, model->A);
1212 cpl_vector_set(vCexp, model_size, model->cexp);
1213 cpl_vector_set(vSigma, model_size, model->sigma);
1214 cpl_vector_set(vB, model_size, model->b);
1215 cpl_vector_set(vC0, model_size, model->c0);
1219 for (
double y = -5; y <= 5; y += 0.1) {
1220 double v = moo_model_eval(model, y);
1222 fprintf(test3,
"%d %f %f\n", x, y, v);
1226 moo_model_delete(model);
1227 cpl_bivector_delete(data);
1231 cpl_vector_set_size(vX, model_size);
1232 cpl_vector_set_size(vA, model_size);
1233 cpl_vector_set_size(vCexp, model_size);
1234 cpl_vector_set_size(vSigma, model_size);
1235 cpl_vector_set_size(vB, model_size);
1236 cpl_vector_set_size(vC0, model_size);
1245 cpl_sprintf(
"%s_%d_model_parameters_step1_medianfilter.txt",
1248 FILE *testm = fopen(filtername,
"w");
1249 fprintf(testm,
"#x A cexp sigma c0 b\n");
1251 for (
int i = 0; i < model_size; i++) {
1252 double x = cpl_vector_get(vX, i);
1253 double A = cpl_vector_get(resA, i);
1254 double cexp = cpl_vector_get(resCexp, i);
1255 double sigma = cpl_vector_get(resSigma, i);
1256 double c0 = cpl_vector_get(resC0, i);
1257 double b = cpl_vector_get(resB, i);
1258 fprintf(testm,
"%f %f %f %f %f %f\n", x, A, cexp, sigma, c0, b);
1261 cpl_free(filtername);
1267 cpl_free(testname1b);
1269 cpl_free(testmodelname);
1270 cpl_free(testmodel3name);
1273 cpl_sprintf(
"%s_%d_localise_step2.txt",
1276 test = fopen(testname,
"w");
1277 fprintf(test,
"#x y flux err qual\n");
1278 cpl_image *rec_fluxes = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
1280 for (
int x = 1; x <= nx; x++) {
1282 double centroid = centroids[(numfib - 1) * nx + x - 1];
1283 double wlo = cpl_image_get(fwlo, x, numfib, &prej);
1284 double wup = cpl_image_get(fwup, x, numfib, &prej);
1285 int yl = floor(centroid - wlo);
1286 int yu = ceil(centroid + wup);
1288 int ysize = yu - yl + 1;
1290 if (!isnan(centroid)) {
1291 cpl_bivector *tmp_data = cpl_bivector_new(ysize);
1292 cpl_vector *tmp_y = cpl_bivector_get_x(tmp_data);
1293 cpl_vector *tmp_flux = cpl_bivector_get_y(tmp_data);
1298 for (
int y = yl; y <= yu; y++) {
1300 hdrl_value flux = hdrl_image_get_pixel(image, x, y, &rej);
1301 fprintf(test,
"%d %d %f %f %d\n", x, y, flux.data,
1304 cpl_vector_set(tmp_y, tmp_size, y);
1305 cpl_vector_set(tmp_flux, tmp_size, flux.data);
1313 double frac = (double)nbrej / (
double)(nbrej + tmp_size);
1315 if (tmp_size > 4 && frac < 0.2) {
1319 cpl_vector_set_size(tmp_flux, tmp_size);
1320 cpl_vector_set_size(tmp_y, tmp_size);
1321 cpl_vector_subtract_scalar(tmp_y, centroid);
1323 int idx = cpl_vector_find(vX, x);
1324 sigma = cpl_vector_get(resSigma, idx);
1325 cexp = cpl_vector_get(resCexp, idx);
1328 _moo_model_fit(tmp_data, &_moo_model_df2, sigma, cexp);
1330 fprintf(testp2,
"%d %f %f %f %f %f %f\n", x, centroid,
1331 model->A, model->sigma, model->c0, model->cexp,
1334 for (
int i = -nstep_lo; i <= nstep_up; i++) {
1335 double y = i * step;
1336 double v = moo_model_eval(model, y);
1337 fprintf(test4,
"%d %f %d %f\n", x, y + centroid, i, v);
1338 cpl_image_set(rec_fluxes, x, i + nstep_lo + 1, v);
1341 moo_model_delete(model);
1344 for (
int i = -nstep_lo; i <= nstep_up; i++) {
1345 double y = i * step;
1346 fprintf(test4,
"%d %f %d %f\n", x, y + centroid, i,
1348 cpl_image_set(rec_fluxes, x, i + nstep_lo + 1, NAN);
1351 cpl_bivector_delete(tmp_data);
1354 fprintf(test4,
"%d %f %d %f\n", x, NAN, 0, NAN);
1355 for (
int i = -nstep_lo; i <= nstep_up; i++) {
1356 cpl_image_set(rec_fluxes, x, i + nstep_lo + 1, NAN);
1361 _moo_median_filter2_model_flat_fibre(rec_fluxes, 10);
1363 char *testmodelp3name =
1364 cpl_sprintf(
"%s_%d_model_step3.txt",
1367 FILE *test5 = fopen(testmodelp3name,
"w");
1368 fprintf(test5,
"#x yr flux\n");
1369 for (
int x = 1; x <= nx; x++) {
1370 for (
int y = 1; y <= ny; y++) {
1372 double v = cpl_image_get(rec_fluxes, x, y, &rej);
1373 fprintf(test5,
"%d %d %f\n", x, y - nstep_lo - 1, v);
1376 cpl_free(testmodelp3name);
1379 cpl_sprintf(
"%s_%d_model_step3.fits",
1382 cpl_image_save(rec_fluxes, testmodelp3name, CPL_TYPE_FLOAT, NULL,
1384 cpl_free(testmodelp3name);
1386 cpl_image_delete(rec_fluxes);
1387 cpl_vector_delete(vX);
1388 cpl_vector_delete(vA);
1389 cpl_vector_delete(vCexp);
1390 cpl_vector_delete(vSigma);
1391 cpl_vector_delete(vC0);
1392 cpl_vector_delete(vB);
1394 cpl_vector_delete(resA);
1395 cpl_vector_delete(resCexp);
1396 cpl_vector_delete(resSigma);
1397 cpl_vector_delete(resB);
1398 cpl_vector_delete(resC0);
1401 cpl_free(testmodel2name);
1407 cpl_free(testmodelp2name);
1419#if (__GNUC__ == 9) && (__GNUC_MINOR__ < 3)
1420#pragma omp parallel shared(nb_fibres, health, nx, ny, image, centroids, \
1421 nstep_lo, nstep_up, step, imgtab, det, \
1422 win_hsize, xstep, fwup, fwlo, extname)
1424#pragma omp parallel default(none) \
1425 shared(nb_fibres, health, nx, ny, image, centroids, nstep_lo, nstep_up, \
1426 step, imgtab, det, win_hsize, xstep, fwup, fwlo, extname)
1431 for (
int i = 1; i <= nb_fibres; i++) {
1432 int h = health[i - 1];
1433 cpl_image *rec_fluxes = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
1435 double *data = cpl_image_get_data(rec_fluxes);
1436 for (
int j = 0; j < nx * ny; j++) {
1450 cpl_msg_info(__func__,
"Modelling profiles %s: fibre %d",
1452 _moo_model_fibre(image, i, centroids, fwlo, fwup, nstep_lo,
1453 nstep_up, step, rec_fluxes, win_hsize, xstep);
1455 _moo_median_filter2_model_flat_fibre(rec_fluxes, 10);
1463 imgtab[i - 1] = rec_fluxes;
1469 cpl_imagelist *cube = cpl_imagelist_new();
1470 for (
int i = 0; i < nb_fibres; i++) {
1471 cpl_imagelist_set(cube, imgtab[i], i);
1476#if MOO_DEBUG_MODEL_FLAT
1478 double crpix2, cd2_2;
1482 cpl_sprintf(
"%s_neg_cube.txt",
1484 FILE *test = fopen(testname,
"w");
1485 fprintf(test,
"#fibnum x y yc flux\n");
1487 for (
int f = 0; f < nb_fibres; f++) {
1488 cpl_image *img = cpl_imagelist_get(cube, f);
1489 int cy = cpl_image_get_size_y(img);
1491 for (
int y = 1; y <= cy; y++) {
1494 for (
int x = 1; x <= nx; x++) {
1495 double centroid = centroids[f * nx + x - 1];
1496 if (!isnan(centroid)) {
1497 double yc = (y - crpix2) * cd2_2 + centroid;
1499 double flux = cpl_image_get(img, x, y, &rej);
1501 fprintf(test,
"%d %d %f %d %f\n", f + 1, x, yc, y,
1513 if (!cpl_errorstate_is_equal(prestate)) {
1514 cpl_msg_error(__func__,
"Error in modelize flat");
1515 cpl_errorstate_dump(prestate, CPL_FALSE, cpl_errorstate_dump_one);
1517 cpl_errorstate_set(prestate);
1548 moo_model_flat_params *params,
1549 const char *filename)
1551 moo_psf *result = NULL;
1553 cpl_errorstate prestate = cpl_errorstate_get();
1555 cpl_ensure(det != NULL, CPL_ERROR_NULL_INPUT, NULL);
1556 cpl_ensure(loc != NULL, CPL_ERROR_NULL_INPUT, NULL);
1557 cpl_ensure(params != NULL, CPL_ERROR_NULL_INPUT, NULL);
1563 cpl_msg_info(__func__,
"Modelizing flat using oversamp %d",
1568 cpl_ensure(fibres_table != NULL, CPL_ERROR_ILLEGAL_INPUT, NULL);
1570 cpl_msg_indent_more();
1571 for (
int i = 1; i <= 2; i++) {
1572 cpl_table_unselect_all(fibres_table);
1574 cpl_table_or_selected_int(fibres_table, MOO_FIBRES_TABLE_SPECTRO,
1575 CPL_EQUAL_TO, fibname);
1576 cpl_table *selected = cpl_table_extract_selected(fibres_table);
1578 cpl_table_get_data_int_const(selected, MOO_FIBRES_TABLE_HEALTH);
1580 for (
int j = 0; j < 3; j++) {
1581 moo_single *det_single =
1584 if (det_single != NULL && loc_single != NULL) {
1585 cpl_msg_info(__func__,
"Modelizing flat for extension %s",
1587 moo_psf_single *psf_single =
1588 _moo_model_flat_single(det_single, loc_single,
1589 result->filename, params, 9, health);
1593 cpl_table_delete(selected);
1595 cpl_msg_indent_less();
1597 if (!cpl_errorstate_is_equal(prestate)) {
1598 cpl_msg_error(__func__,
"Error in model flat");
1599 cpl_errorstate_dump(prestate, CPL_FALSE, cpl_errorstate_dump_one);
1601 cpl_errorstate_set(prestate);
#define MOO_BADPIX_OUTSIDE_DATA_RANGE
#define MOO_BADPIX_COSMETIC
moo_single * moo_det_load_single(moo_det *self, moo_detector_type type, int num, int level)
Load the type part in DET and return it.
const char * moo_detector_get_extname(moo_detector_type type, int ntas)
Get the extension name of a detector.
cpl_image * moo_loc_single_get_f_wup(moo_loc_single *self)
Get image of width low.
cpl_image * moo_loc_single_get_f_centroids(moo_loc_single *self)
Get image of fit centroids.
cpl_image * moo_loc_single_get_f_wlo(moo_loc_single *self)
Get image of width low.
moo_loc_single * moo_loc_get_single(moo_loc *self, moo_detector_type type, int ntas)
Get the type part in LOC and return it.
cpl_table * moo_loc_get_fibre_table(moo_loc *self)
Get the FIBRE TABLE in LOC.
cpl_error_code moo_psf_single_set_cube_ref(moo_psf_single *self, double crpix2, double cd2_2)
Set cube parameters.
cpl_error_code moo_psf_single_get_cube_ref(moo_psf_single *self, double *crpix2, double *cd2_2)
Set cube parameters.
moo_psf_single * moo_psf_single_create(const char *filename, const char *extname)
Create a new moo_psf_single from the given PSF filename.
cpl_error_code moo_psf_add_single(moo_psf *self, moo_psf_single *single, moo_detector_type type, int ntas)
Add PSF_SINGLE extension to PSF filename and update moo_psf structure.
moo_psf * moo_psf_create(const char *filename)
Create a new empty PSF filename.
hdrl_image * moo_single_get_image(moo_single *self)
Get the IMAGE part (DATA,ERR) of single DET.
moo_psf * moo_model_flat(moo_det *det, moo_loc *loc, moo_model_flat_params *params, const char *filename)
To extract 2D FF and model it, resulting MASTER_FLAT.
cpl_vector * moo_savgol_filter(cpl_vector *v, int window_length, int polyorder)
Apply a Savitzky-Golay filter to a vector.
double moo_vector_get_min(const cpl_vector *v, int *flags)
Find minimum values in a vector using flags.
double moo_vector_get_max(const cpl_vector *v, int *flags)
Find maximum values in a vector using flags.
cpl_error_code moo_tchebychev_fit(cpl_bivector *data, int *flag, int degree, double xmin, double xmax, double ymin, double ymax)
Computes Tchebitchev transformation of data.
cpl_vector * moo_median_filter(cpl_vector *v, int winhsize)
Apply a median filter to a vector.
double moo_vector_get_percentile(cpl_vector *v, double f)
Get percentile of input vector.