MOONS Pipeline Reference Manual 0.13.1
moo_model_flat.c
1/*
2 * This file is part of the MOONS Pipeline
3 * Copyright (C) 2002-2016 European Southern Observatory
4 *
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
18 */
19
20#ifdef HAVE_CONFIG_H
21#include <config.h>
22#endif
23
24/*-----------------------------------------------------------------------------
25 Includes
26 -----------------------------------------------------------------------------*/
27#include <math.h>
28#include <string.h>
29#include <cpl.h>
30#include <hdrl.h>
31#include <gsl/gsl_vector.h>
32#include <gsl/gsl_blas.h>
33#include <gsl/gsl_multifit_nlin.h>
34#include "moo_utils.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"
41#ifdef _OPENMP
42#include <omp.h>
43#endif
44/*----------------------------------------------------------------------------*/
49/*----------------------------------------------------------------------------*/
52/*-----------------------------------------------------------------------------
53 Function codes
54 -----------------------------------------------------------------------------*/
55double
56moo_image_get_linear_y_interpolated(cpl_image *image, double y, int x, int *rej)
57{
58 double flux = NAN;
59
60 cpl_ensure(image != NULL, CPL_ERROR_NULL_INPUT, NAN);
61 cpl_ensure(rej != NULL, CPL_ERROR_NULL_INPUT, NAN);
62
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);
67
68 int y1 = (int)floor(y);
69 int y2 = y1 + 1;
70
71 double flux1 = cpl_image_get(image, x, y1, rej);
72
73 double flux2 = cpl_image_get(image, x, y2, rej);
74
75 if (*rej == CPL_BINARY_0) {
76 flux = flux1 + (y - y1) * (flux2 - flux1);
77 }
78 return flux;
79}
80
81double
82moo_image_get_pixel_y_interpolated(cpl_image *image,
83 double y,
84 int x,
85 double step,
86 int *rej)
87{
88 double flux = NAN;
89
90 cpl_ensure(image != NULL, CPL_ERROR_NULL_INPUT, NAN);
91 cpl_ensure(rej != NULL, CPL_ERROR_NULL_INPUT, NAN);
92
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);
97
98 double ystart = y - step / 2.0;
99 double ystop = y + step / 2.0;
100
101 int pix1 = (int)floor(ystart + 0.5);
102 int pix2 = (int)floor(ystop + 0.5);
103
104 if (pix1 == pix2) {
105 flux = cpl_image_get(image, x, pix1, rej);
106 }
107 else {
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);
112
113 flux = flux1 * frac1 + flux2 * frac2;
114 // cpl_msg_info("test","get %f * %d (%f) + %f * %d (%f) = %f",frac1,pix1,flux1,frac2,pix2,flux2,flux);
115 }
116
117 if (*rej == CPL_BINARY_1) {
118 flux = NAN;
119 }
120
121 return flux;
122}
123
124double
125moo_image_get_pixel2_y_interpolated(cpl_image *image, double y, int x, int *rej)
126{
127 double flux = NAN;
128
129 cpl_ensure(image != NULL, CPL_ERROR_NULL_INPUT, NAN);
130 cpl_ensure(rej != NULL, CPL_ERROR_NULL_INPUT, NAN);
131
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);
136
137 int pix = (int)floor(y + 0.5);
138
139 flux = cpl_image_get(image, x, pix, rej);
140
141 if (*rej == CPL_BINARY_1) {
142 flux = NAN;
143 }
144
145 return flux;
146}
147
148double
149moo_image_get_pixel_yerror_interpolated(cpl_image *error, double y, int x)
150{
151 double err = NAN;
152
153 cpl_ensure(error != NULL, CPL_ERROR_NULL_INPUT, NAN);
154
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);
159
160 int y1 = (int)floor(y + 0.5);
161
162 int rej;
163 err = cpl_image_get(error, x, y1, &rej);
164
165 return err;
166}
167
168double
169moo_image_get_linear_yerror_interpolated(cpl_image *error, double y, int x)
170{
171 double err = NAN;
172
173 cpl_ensure(error != NULL, CPL_ERROR_NULL_INPUT, NAN);
174
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);
179
180 int y1 = (int)floor(y);
181 int y2 = y1 + 1;
182
183 int rej;
184 double err1 = cpl_image_get(error, x, y1, &rej);
185
186 double err2 = cpl_image_get(error, x, y2, &rej);
187
188 err = sqrt((y - y1) * err1 * (y - y1) * err1 +
189 (y2 - y) * err2 * (y2 - y) * err2);
190
191 return err;
192}
193
194static void
195_moo_fit_model_flat_fibre(cpl_image *rec_fluxes, int deg_poly)
196{
197 int nx = cpl_image_get_size_x(rec_fluxes);
198 int ny = cpl_image_get_size_y(rec_fluxes);
199
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);
203
204 int *flag = cpl_calloc(nx, sizeof(int));
205
206 for (int x = 0; x < nx; x++) {
207 cpl_vector_set(xpos, x, x + 1);
208 }
209
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);
215 if (isnan(val)) {
216 flag[x] = 1;
217 }
218 else {
219 flag[x] = 0;
220 }
221 }
222 double min = moo_vector_get_min(fluxes, flag);
223 double max = moo_vector_get_max(fluxes, flag);
224
225 moo_tchebychev_fit(points, flag, deg_poly, 1, nx, min, max);
226
227 for (int x = 0; x < nx; x++) {
228 double val = cpl_vector_get(fluxes, x);
229 data[x + y * nx] = val;
230 }
231 }
232
233 cpl_bivector_delete(points);
234 cpl_free(flag);
235}
236
237static void
238_moo_smooth_model_flat_fibre(cpl_image *rec_fluxes, double s)
239{
240 int nx = cpl_image_get_size_x(rec_fluxes);
241 int ny = cpl_image_get_size_y(rec_fluxes);
242
243 cpl_vector *fluxes = cpl_vector_new(nx);
244
245 double *data = cpl_image_get_data(rec_fluxes);
246 for (int y = 0; y < ny; y++) {
247 int nbnan = 0;
248 for (int x = 0; x < nx; x++) {
249 double val = data[x + y * nx];
250 if (isnan(val)) {
251 nbnan++;
252 }
253 cpl_vector_set(fluxes, x, val);
254 }
255
256 cpl_vector *filter = cpl_vector_new(nx - nbnan);
257
258 int filteridx = 0;
259 for (int x = 0; x < nx; x++) {
260 double val = cpl_vector_get(fluxes, x);
261
262 if (!isnan(val)) {
263 cpl_vector_set(filter, filteridx, val);
264 filteridx++;
265 }
266 }
267
268 cpl_vector *res = moo_hpfilter(filter, s);
269
270 nbnan = 0;
271
272 for (int x = 0; x < nx; x++) {
273 double val = cpl_vector_get(fluxes, x);
274 if (!isnan(val)) {
275 val = cpl_vector_get(res, x - nbnan);
276 }
277 else {
278 nbnan++;
279 }
280 data[x + y * nx] = val;
281 }
282 cpl_vector_delete(filter);
283 cpl_vector_delete(res);
284 }
285
286 cpl_vector_delete(fluxes);
287}
288
289static void
290_moo_median_filter_model_flat_fibre(cpl_image *rec_fluxes,
291 int winhsize,
292 int kappa)
293{
294 int nx = cpl_image_get_size_x(rec_fluxes);
295 int ny = cpl_image_get_size_y(rec_fluxes);
296
297 cpl_mask *mask = cpl_mask_new(nx, ny);
298
299 cpl_mask *orig = cpl_image_get_bpm(rec_fluxes);
300
301 cpl_vector *fluxes = cpl_vector_new(nx);
302 int window_size = winhsize * 2 + 1;
303
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++) {
307 int rej;
308 double val = cpl_image_get(rec_fluxes, x, y, &rej);
309 cpl_vector_set(fluxes, x - 1, val);
310 }
311
312 for (int x = 0; x < window_size; x++) {
313 double val = cpl_vector_get(fluxes, x);
314 cpl_vector_set(winflux, x, val);
315 }
316 double median = moo_vector_get_percentile(winflux, 0.5);
317 double p1 = moo_vector_get_percentile(winflux, 0.1585);
318 double p2 = moo_vector_get_percentile(winflux, 0.8415);
319 double sigma = (p2 - p1) / 2;
320
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);
326 }
327 }
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);
332 }
333 double val = cpl_vector_get(winflux, winhsize);
334 if (!isnan(val)) {
335 median = moo_vector_get_percentile(winflux, 0.5);
336 p1 = moo_vector_get_percentile(winflux, 0.1585);
337 p2 = moo_vector_get_percentile(winflux, 0.8415);
338 sigma = (p2 - p1) / 2;
339
340 double diff = val - median;
341 if (fabs(diff) > kappa * sigma) {
342 cpl_mask_set(mask, x + 1, y, CPL_BINARY_1);
343 }
344 }
345 }
346 }
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);
353}
354
355static void
356_moo_median_filter2_model_flat_fibre(cpl_image *rec_fluxes, int winhsize)
357{
358 int nx = cpl_image_get_size_x(rec_fluxes);
359 int ny = cpl_image_get_size_y(rec_fluxes);
360
361 cpl_vector *fluxes = cpl_vector_new(nx);
362
363 for (int y = 1; y <= ny; y++) {
364 for (int x = 1; x <= nx; x++) {
365 int rej;
366 double val = cpl_image_get(rec_fluxes, x, y, &rej);
367 cpl_vector_set(fluxes, x - 1, val);
368 }
369
370 cpl_vector *res = moo_median_filter(fluxes, winhsize);
371
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);
375 }
376
377 cpl_vector_delete(res);
378 }
379 cpl_vector_delete(fluxes);
380}
381
382static void
383_moo_savgol_model_flat_fibre(cpl_image *rec_fluxes,
384 int polyorder,
385 int window_length)
386{
387 int nx = cpl_image_get_size_x(rec_fluxes);
388 int ny = cpl_image_get_size_y(rec_fluxes);
389
390 cpl_vector *fluxes = cpl_vector_new(nx);
391
392 for (int y = 1; y <= ny; y++) {
393 for (int x = 1; x <= nx; x++) {
394 int rej;
395 double val = cpl_image_get(rec_fluxes, x, y, &rej);
396
397 if (rej != 0) {
398 val = NAN;
399 }
400 cpl_vector_set(fluxes, x - 1, val);
401 }
402 if (y == 1) {
403 char *testname = cpl_sprintf("befsavgol.txt");
404 FILE *test = fopen(testname, "w");
405 fprintf(test, "#x f\n");
406
407 for (int x = 0; x < nx; x++) {
408 double val = cpl_vector_get(fluxes, x);
409 fprintf(test, "%d %f\n", x, val);
410 }
411
412 fclose(test);
413 cpl_free(testname);
414 }
415 cpl_vector *res = moo_savgol_filter(fluxes, window_length, polyorder);
416 if (y == 1) {
417 char *testname = cpl_sprintf("savgol.txt");
418 FILE *test = fopen(testname, "w");
419 fprintf(test, "#x f\n");
420
421 for (int x = 0; x < nx; x++) {
422 double val = cpl_vector_get(res, x);
423 fprintf(test, "%d %f\n", x, val);
424 }
425
426 fclose(test);
427 cpl_free(testname);
428 }
429
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);
433 }
434 cpl_vector_delete(res);
435 }
436
437 cpl_vector_delete(fluxes);
438}
439
440static void
441_moo_extract_error_model_flat_fibre(hdrl_image *image,
442 int numfib,
443 double *centroids,
444 int nstep_lo,
445 int nstep_up,
446 double step,
447 cpl_image *rec_fluxes)
448{
449 int nx = hdrl_image_get_size_x(image);
450 cpl_image *error = hdrl_image_get_error(image);
451
452 for (int x = 1; x <= nx; x++) {
453 double centroid = centroids[(numfib - 1) * nx + x - 1];
454
455 if (!isnan(centroid)) {
456 for (int y = -nstep_lo; y <= nstep_up; y++) {
457 double yc = centroid + y * step;
458 double flux =
459 moo_image_get_pixel_yerror_interpolated(error, yc, x);
460 cpl_image_set(rec_fluxes, x, y + nstep_lo + 1, flux);
461 }
462 }
463 else {
464 for (int y = -nstep_lo; y <= nstep_up; y++) {
465 cpl_image_set(rec_fluxes, x, y + nstep_lo + 1, NAN);
466 }
467 }
468 }
469}
470
471static void
472_moo_rectify_fibre(hdrl_image *image,
473 int numfib,
474 double *centroids,
475 int nstep_lo,
476 int nstep_up,
477 double step,
478 cpl_image *rec_fluxes)
479{
480 int nx = hdrl_image_get_size_x(image);
481 cpl_image *data = hdrl_image_get_image(image);
482
483 for (int x = 1; x <= nx; x++) {
484 double centroid = centroids[(numfib - 1) * nx + x - 1];
485
486 if (!isnan(centroid)) {
487 for (int y = -nstep_lo; y <= nstep_up; y++) {
488 double yc = centroid + y * step;
489 int rej = 0;
490 // double flux = moo_image_get_linear_y_interpolated( data, yc, x, &rej);
491 // double flux = moo_image_get_pixel_y_interpolated( data, yc, x, step, &rej);
492 double flux =
493 moo_image_get_pixel2_y_interpolated(data, yc, x, &rej);
494 cpl_image_set(rec_fluxes, x, y + nstep_lo + 1, flux);
495 }
496 }
497 else {
498 for (int y = -nstep_lo; y <= nstep_up; y++) {
499 cpl_image_set(rec_fluxes, x, y + nstep_lo + 1, NAN);
500 }
501 }
502 }
503}
504
505static int
506_moo_model_f(const gsl_vector *px, void *data, gsl_vector *f)
507{
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);
512
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);
518
519 moo_model *m = moo_model_new(sigma, cexp, A, c0, b);
520
521 size_t i;
522
523 for (i = 0; i < n; i++) {
524 /* Model Yi = A * exp(-lambda * i) + b */
525 double x = cpl_vector_get(vx, i);
526 double y = cpl_vector_get(vy, i);
527
528 double Yx = moo_model_eval(m, x);
529
530 gsl_vector_set(f, i, Yx - y);
531 }
532 moo_model_delete(m);
533
534 return GSL_SUCCESS;
535}
536
537static int
538_moo_model_df(const gsl_vector *px, void *data, gsl_matrix *J)
539{
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);
543
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);
549
550 moo_model *m = moo_model_new(sigma, cexp, A, c0, b);
551 size_t i;
552
553 for (i = 0; i < n; i++) {
554 /* Jacobian matrix J(i,j) = dfi / dxj, */
555 /* where fi = (Yi - yi)/sigma[i], */
556 /* Yi = A * exp(-lambda * i) + b */
557 /* and the xj are the parameters (A,lambda,b) */
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);
563
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);
569 }
570 moo_model_delete(m);
571
572 return GSL_SUCCESS;
573}
574
575static int
576_moo_model_df2(const gsl_vector *px, void *data, gsl_matrix *J)
577{
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);
581
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);
587
588 moo_model *m = moo_model_new(sigma, cexp, A, c0, b);
589 size_t i;
590
591 for (i = 0; i < n; i++) {
592 /* Jacobian matrix J(i,j) = dfi / dxj, */
593 /* where fi = (Yi - yi)/sigma[i], */
594 /* Yi = A * exp(-lambda * i) + b */
595 /* and the xj are the parameters (A,lambda,b) */
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);
599
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);
605 }
606 moo_model_delete(m);
607
608 return GSL_SUCCESS;
609}
610
611static moo_model *
612_moo_model_fit(cpl_bivector *data,
613 int (*fit_func)(const gsl_vector *, void *, gsl_matrix *),
614 double isigma,
615 double icexp)
616{
617 moo_model *model = NULL;
618
619 if (data != 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;
625 double b = min;
626 const gsl_multifit_fdfsolver_type *T = gsl_multifit_fdfsolver_lmsder;
627 gsl_multifit_fdfsolver *s;
628 int info;
629
630 const size_t p = 5;
631
632 gsl_matrix *J = gsl_matrix_alloc(N, p);
633 gsl_matrix *covar = gsl_matrix_alloc(p, p);
634
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);
638
639 const double xtol = 1e-9;
640 const double gtol = 1e-9;
641 const double ftol = 0.0;
642
643 f.f = &_moo_model_f;
644 f.df = fit_func;
645 f.n = N;
646 f.p = p;
647 f.params = data;
648
649 s = gsl_multifit_fdfsolver_alloc(T, N, p);
650
651 /* initialize solver with starting point */
652 gsl_multifit_fdfsolver_set(s, &f, &x.vector);
653
654 /* compute initial residual norm */
655 gsl_multifit_fdfsolver_residual(s);
656
657 /* solve the system with a maximum of 50 iterations */
658 gsl_multifit_fdfsolver_driver(s, 50, xtol, gtol, ftol, &info);
659
660 gsl_multifit_fdfsolver_jac(s, J);
661 gsl_multifit_covar(J, 0.0, covar);
662
663 /*
664 double chi = gsl_blas_dnrm2(res_f);
665
666 #define FIT(i) gsl_vector_get(s->x, i)
667 #define ERR(i) sqrt(gsl_matrix_get(covar,i,i))
668 {
669 double dof = N - p;
670 double c = GSL_MAX_DBL(1, chi / sqrt(dof));
671
672 fprintf(stderr, "chisq/dof = %g\n", pow(chi, 2.0) / dof);
673
674 fprintf (stderr, "sigma = %.5f +/- %.5f\n", FIT(0), c*ERR(0));
675 fprintf (stderr, "cexp = %.5f +/- %.5f\n", FIT(1), c*ERR(1));
676 fprintf (stderr, "A = %.5f +/- %.5f\n", FIT(2), c*ERR(2));
677 fprintf (stderr, "c0 = %.5f +/- %.5f\n", FIT(3), c*ERR(3));
678 fprintf (stderr, "B = %.5f +/- %.5f\n", FIT(4), c*ERR(4));
679 }
680
681 fprintf (stderr, "status = %s\n", gsl_strerror (status));
682 */
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);
688
689 model = moo_model_new(sigma, cexp, A, c0, b);
690 gsl_multifit_fdfsolver_free(s);
691 gsl_matrix_free(covar);
692 gsl_matrix_free(J);
693 }
694 return model;
695}
696
697moo_model *
698moo_model_new(double sigma, double cexp, double A, double c0, double b)
699{
700 moo_model *res = cpl_calloc(1, sizeof(moo_model));
701
702 res->sigma = sigma;
703 res->cexp = cexp;
704 res->A = A;
705 res->c0 = c0;
706 res->b = b;
707
708 return res;
709}
710
711void
712moo_model_delete(moo_model *self)
713{
714 if (self != NULL) {
715 cpl_free(self);
716 }
717}
718
719double
720moo_model_eval(moo_model *model, double x)
721{
722 double res = NAN;
723
724 cpl_ensure(model != NULL, CPL_ERROR_NULL_INPUT, res);
725
726 double sigma = model->sigma;
727 double cexp = model->cexp;
728 double A = model->A;
729 double c0 = model->c0;
730 double b = model->b;
731
732 res = (A / (sqrt(CPL_MATH_2PI) * sigma)) *
733 exp(-0.5 * pow((fabs(x - c0) / sigma), cexp)) +
734 b;
735
736 return res;
737}
738
739double
740moo_model_dA_eval(moo_model *model, double x)
741{
742 double res = NAN;
743
744 cpl_ensure(model != NULL, CPL_ERROR_NULL_INPUT, res);
745
746 double sigma = model->sigma;
747 double cexp = model->cexp;
748 double c0 = model->c0;
749
750 double absxc0 = fabs(x - c0);
751
752 res = (1 / (sqrt(CPL_MATH_2PI) * sigma)) *
753 exp(-0.5 * pow((absxc0 / sigma), cexp));
754
755 return res;
756}
757
758double
759moo_model_dsigma_eval(moo_model *model, double x)
760{
761 double res = NAN;
762
763 cpl_ensure(model != NULL, CPL_ERROR_NULL_INPUT, res);
764
765 double sigma = model->sigma;
766 double cexp = model->cexp;
767 double c0 = model->c0;
768 double A = model->A;
769
770 double absxc0 = fabs(x - c0);
771
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));
776 return res;
777}
778
779double
780moo_model_dcexp_eval(moo_model *model, double x)
781{
782 double res = NAN;
783
784 cpl_ensure(model != NULL, CPL_ERROR_NULL_INPUT, res);
785
786 double sigma = model->sigma;
787 double cexp = model->cexp;
788 double c0 = model->c0;
789 double A = model->A;
790
791 if (x == c0) {
792 res = 0;
793 }
794 else {
795 double absxc0 = fabs(x - c0);
796
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));
801 }
802 return res;
803}
804
805double
806moo_model_dc0_eval(moo_model *model, double x)
807{
808 double res = NAN;
809
810 cpl_ensure(model != NULL, CPL_ERROR_NULL_INPUT, res);
811
812 double sigma = model->sigma;
813 double cexp = model->cexp;
814 double c0 = model->c0;
815 double A = model->A;
816
817 double absxc0 = fabs(x - c0);
818
819 if (x == c0) {
820 res = 0;
821 }
822 else {
823 res = (A * cexp * pow(sigma, -cexp - 1) *
824 exp(-pow(absxc0, cexp) / (2 * pow(sigma, cexp))) *
825 pow(absxc0, cexp)) /
826 (2 * sqrt(CPL_MATH_2PI) * (x - c0));
827 }
828 return res;
829}
830
831static void
832_moo_model_fibre(hdrl_image *image,
833 int numfib,
834 double *centroids,
835 cpl_image *fwlo,
836 cpl_image *fwup,
837 int nstep_lo,
838 int nstep_up,
839 double step,
840 cpl_image *rec_fluxes,
841 int win_hsize,
842 int xstep)
843{
844 int nx = hdrl_image_get_size_x(image);
845
846 cpl_vector *vX = cpl_vector_new(nx);
847 cpl_vector *vCexp = cpl_vector_new(nx);
848 cpl_vector *vSigma = cpl_vector_new(nx);
849
850 int model_size = 0;
851
852 for (int x = win_hsize + 1; x <= nx - win_hsize; x += xstep) {
853 cpl_bivector *data = NULL;
854
855 for (int x2 = x - win_hsize; x2 <= x + win_hsize; x2++) {
856 int prej;
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);
862
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);
867
868 int tmp_size = 0;
869 int nbrej = 0;
870 if (!isnan(centroid)) {
871 for (int y = yl; y <= yu; y++) {
872 int rej = 0;
873 hdrl_value flux = hdrl_image_get_pixel(image, x2, y, &rej);
874 if (rej == 0) {
875 cpl_vector_set(tmp_y, tmp_size, y);
876 cpl_vector_set(tmp_flux, tmp_size, flux.data);
877 tmp_size++;
878 }
879 else {
880 nbrej++;
881 }
882 }
883 }
884
885 double frac = 1 - (double)nbrej / (double)(tmp_size + nbrej);
886
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);
892
893 cpl_vector_subtract_scalar(tmp_flux, min);
894
895 double sum = cpl_vector_get_sum(tmp_flux);
896
897 if (sum > 0) {
898 cpl_vector_divide_scalar(tmp_flux, sum);
899
900 if (data == NULL) {
901 data = cpl_bivector_duplicate(tmp_data);
902 }
903 else {
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));
914 }
915 }
916 }
917 }
918 cpl_bivector_delete(tmp_data);
919 }
920
921 moo_model *model = _moo_model_fit(data, &_moo_model_df, 1., 2.);
922
923 if (model != NULL) {
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);
927 model_size++;
928 }
929
930 moo_model_delete(model);
931 cpl_bivector_delete(data);
932 }
933
934 cpl_vector_set_size(vX, model_size);
935 cpl_vector_set_size(vCexp, model_size);
936 cpl_vector_set_size(vSigma, model_size);
937
938 cpl_vector *resCexp = moo_median_filter(vCexp, 5);
939 cpl_vector *resSigma = moo_median_filter(vSigma, 5);
940
941 for (int x = 1; x <= nx; x++) {
942 int prej;
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);
948
949 int ysize = yu - yl + 1;
950
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);
955
956 int tmp_size = 0;
957 int nbrej = 0;
958
959 for (int y = yl; y <= yu; y++) {
960 int rej = 0;
961 hdrl_value flux = hdrl_image_get_pixel(image, x, y, &rej);
962 if (rej == 0) {
963 cpl_vector_set(tmp_y, tmp_size, y);
964 cpl_vector_set(tmp_flux, tmp_size, flux.data);
965 tmp_size++;
966 }
967 else {
968 nbrej++;
969 }
970 }
971 double frac = (double)nbrej / (double)(nbrej + tmp_size);
972
973 if (tmp_size > 4 && frac < 0.2) {
974 double sigma = 1.0;
975 double cexp = 2.0;
976
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);
980
981 int idx = cpl_vector_find(vX, x);
982 sigma = cpl_vector_get(resSigma, idx);
983 cexp = cpl_vector_get(resCexp, idx);
984
985 moo_model *model =
986 _moo_model_fit(tmp_data, &_moo_model_df2, sigma, cexp);
987 for (int i = -nstep_lo; i <= nstep_up; i++) {
988 double y = i * step;
989 double v = moo_model_eval(model, y);
990
991 cpl_image_set(rec_fluxes, x, i + nstep_lo + 1, v);
992 }
993
994 moo_model_delete(model);
995 }
996 else {
997 for (int i = -nstep_lo; i <= nstep_up; i++) {
998 cpl_image_set(rec_fluxes, x, i + nstep_lo + 1, NAN);
999 }
1000 }
1001 cpl_bivector_delete(tmp_data);
1002 }
1003 else {
1004 for (int i = -nstep_lo; i <= nstep_up; i++) {
1005 cpl_image_set(rec_fluxes, x, i + nstep_lo + 1, NAN);
1006 }
1007 }
1008 }
1009
1010 cpl_vector_delete(vX);
1011 cpl_vector_delete(vCexp);
1012 cpl_vector_delete(vSigma);
1013
1014 cpl_vector_delete(resCexp);
1015 cpl_vector_delete(resSigma);
1016}
1017
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,
1023 int deg_poly,
1024 const int *health)
1025{
1026 cpl_errorstate prestate = cpl_errorstate_get();
1027
1028 cpl_ensure(filename != NULL, CPL_ERROR_NULL_INPUT, NULL);
1029 cpl_ensure(det != NULL, CPL_ERROR_NULL_INPUT, NULL);
1030
1031 cpl_image *fcentroids = moo_loc_single_get_f_centroids(loc);
1032 cpl_image *fwlo = moo_loc_single_get_f_wlo(loc);
1033 cpl_image *fwup = moo_loc_single_get_f_wup(loc);
1034 hdrl_image *image = moo_single_get_image(det);
1035
1036 int nb_fibres = cpl_image_get_size_y(fcentroids);
1037 int nx = cpl_image_get_size_x(fcentroids);
1038
1039 cpl_image_reject_value(fwlo, CPL_VALUE_NAN);
1040 cpl_image_reject_value(fwup, CPL_VALUE_NAN);
1041
1042 double max_fwlo = cpl_image_get_max(fwlo);
1043 double max_fwup = cpl_image_get_max(fwup);
1044
1045 cpl_msg_indent_more();
1046 cpl_msg_info("moo_model_flat", "Use Y localisation limits: %.3f,%.3f pix",
1047 -max_fwlo, max_fwup);
1048
1049 int nstep_lo = (int)ceil(max_fwlo * params->oversamp);
1050 int nstep_up = (int)ceil(max_fwup * params->oversamp);
1051
1052 int ny = nstep_lo + nstep_up + 1;
1053
1054 moo_psf_single *res = moo_psf_single_create(filename, det->extname);
1055
1056 double *centroids = cpl_image_get_data(fcentroids);
1057 double step = 1.0 / params->oversamp;
1058
1059 cpl_image **imgtab = cpl_calloc(nb_fibres, sizeof(cpl_image *));
1060 const char *extname = moo_detector_get_extname(det->type, det->ntas);
1061
1062 int win_hsize = params->winhsize;
1063 int xstep = params->xstep;
1064
1065 cpl_msg_info("moo_model_flat",
1066 "Fitting profile shape parameters using step %d and window "
1067 "half size %d along X axis",
1068 xstep, win_hsize);
1069 cpl_msg_indent_less();
1070
1071#if MOO_DEBUG_MODEL_FLAT
1072 int test_fibers[] = { 10, 250, 500 };
1073 for (int ti = 306; ti <= 308; ti++) {
1074 int numfib = ti;
1075
1076 char *testname =
1077 cpl_sprintf("%s_%d_localise_step1.txt",
1078 moo_detector_get_extname(det->type, det->ntas), numfib);
1079 char *testname1b =
1080 cpl_sprintf("%s_%d_localise_align_centroid_step1.txt",
1081 moo_detector_get_extname(det->type, det->ntas), numfib);
1082 char *testmodelname =
1083 cpl_sprintf("%s_%d_model_parameters_step1.txt",
1084 moo_detector_get_extname(det->type, det->ntas), numfib);
1085 char *testmodel2name =
1086 cpl_sprintf("%s_%d_model_step1.txt",
1087 moo_detector_get_extname(det->type, det->ntas), numfib);
1088 char *testmodel3name =
1089 cpl_sprintf("%s_%d_model_step2.txt",
1090 moo_detector_get_extname(det->type, det->ntas), numfib);
1091 char *testmodelp2name =
1092 cpl_sprintf("%s_%d_model_parameters_step2.txt",
1093 moo_detector_get_extname(det->type, det->ntas), numfib);
1094 FILE *test = fopen(testname, "w");
1095 fprintf(test, "#x x2 y flux err qual\n");
1096
1097 FILE *test1 = fopen(testname1b, "w");
1098 fprintf(test1, "#x x2 y flux centroid min sum\n");
1099
1100 FILE *test2 = fopen(testmodelname, "w");
1101 fprintf(test2, "#x centroid A sigma c0 cexp b\n");
1102
1103 FILE *test3 = fopen(testmodel2name, "w");
1104 fprintf(test3, "#x y flux\n");
1105
1106 FILE *test4 = fopen(testmodel3name, "w");
1107 fprintf(test4, "#x y yr flux\n");
1108
1109 FILE *testp2 = fopen(testmodelp2name, "w");
1110 fprintf(testp2, "#x centroid A sigma c0 cexp b\n");
1111
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);
1118
1119 int model_size = 0;
1120
1121 for (int x = win_hsize + 1; x <= nx - win_hsize; x += xstep) {
1122 cpl_bivector *data = NULL;
1123
1124 for (int x2 = x - win_hsize; x2 <= x + win_hsize; x2++) {
1125 int prej;
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);
1131
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);
1136
1137 int tmp_size = 0;
1138 int nbrej = 0;
1139
1140 if (!isnan(centroid)) {
1141 for (int y = yl; y <= yu; y++) {
1142 int rej = 0;
1143 hdrl_value flux =
1144 hdrl_image_get_pixel(image, x2, y, &rej);
1145 if (rej == 0) {
1146 fprintf(test, "%d %d %d %f %f %d\n", x, x2, y,
1147 flux.data, flux.error, rej);
1148
1149 cpl_vector_set(tmp_y, tmp_size, y);
1150 cpl_vector_set(tmp_flux, tmp_size, flux.data);
1151 tmp_size++;
1152 }
1153 else {
1154 nbrej++;
1155 }
1156 }
1157 }
1158 double frac = 1 - (double)nbrej / (double)(tmp_size + nbrej);
1159
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);
1163
1164 double min = cpl_vector_get_min(tmp_flux);
1165
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);
1169
1170 if (sum > 0) {
1171 cpl_vector_divide_scalar(tmp_flux, sum);
1172
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);
1178 }
1179
1180 if (data == NULL) {
1181 data = cpl_bivector_duplicate(tmp_data);
1182 }
1183 else {
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,
1194 i - size));
1195 }
1196 }
1197 }
1198 }
1199
1200 cpl_bivector_delete(tmp_data);
1201 }
1202
1203 if (data != NULL) {
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,
1209 model->b);
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);
1216
1217 model_size++;
1218
1219 for (double y = -5; y <= 5; y += 0.1) {
1220 double v = moo_model_eval(model, y);
1221
1222 fprintf(test3, "%d %f %f\n", x, y, v);
1223 }
1224 }
1225
1226 moo_model_delete(model);
1227 cpl_bivector_delete(data);
1228 }
1229 }
1230
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);
1237
1238 cpl_vector *resA = moo_median_filter(vA, 5);
1239 cpl_vector *resCexp = moo_median_filter(vCexp, 5);
1240 cpl_vector *resSigma = moo_median_filter(vSigma, 5);
1241 cpl_vector *resB = moo_median_filter(vB, 5);
1242 cpl_vector *resC0 = moo_median_filter(vC0, 5);
1243 {
1244 char *filtername =
1245 cpl_sprintf("%s_%d_model_parameters_step1_medianfilter.txt",
1246 moo_detector_get_extname(det->type, det->ntas),
1247 numfib);
1248 FILE *testm = fopen(filtername, "w");
1249 fprintf(testm, "#x A cexp sigma c0 b\n");
1250
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);
1259 }
1260 fclose(testm);
1261 cpl_free(filtername);
1262 }
1263
1264 fclose(test);
1265 cpl_free(testname);
1266 fclose(test1);
1267 cpl_free(testname1b);
1268 fclose(test2);
1269 cpl_free(testmodelname);
1270 cpl_free(testmodel3name);
1271
1272 testname =
1273 cpl_sprintf("%s_%d_localise_step2.txt",
1274 moo_detector_get_extname(det->type, det->ntas), numfib);
1275
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);
1279
1280 for (int x = 1; x <= nx; x++) {
1281 int prej;
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);
1287
1288 int ysize = yu - yl + 1;
1289
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);
1294
1295 int tmp_size = 0;
1296 int nbrej = 0;
1297
1298 for (int y = yl; y <= yu; y++) {
1299 int rej = 0;
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,
1302 flux.error, rej);
1303 if (rej == 0) {
1304 cpl_vector_set(tmp_y, tmp_size, y);
1305 cpl_vector_set(tmp_flux, tmp_size, flux.data);
1306 tmp_size++;
1307 }
1308 else {
1309 nbrej++;
1310 }
1311 }
1312
1313 double frac = (double)nbrej / (double)(nbrej + tmp_size);
1314
1315 if (tmp_size > 4 && frac < 0.2) {
1316 double sigma = 1.0;
1317 double cexp = 2.0;
1318
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);
1322
1323 int idx = cpl_vector_find(vX, x);
1324 sigma = cpl_vector_get(resSigma, idx);
1325 cexp = cpl_vector_get(resCexp, idx);
1326
1327 moo_model *model =
1328 _moo_model_fit(tmp_data, &_moo_model_df2, sigma, cexp);
1329
1330 fprintf(testp2, "%d %f %f %f %f %f %f\n", x, centroid,
1331 model->A, model->sigma, model->c0, model->cexp,
1332 model->b);
1333
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);
1339 }
1340
1341 moo_model_delete(model);
1342 }
1343 else {
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,
1347 NAN);
1348 cpl_image_set(rec_fluxes, x, i + nstep_lo + 1, NAN);
1349 }
1350 }
1351 cpl_bivector_delete(tmp_data);
1352 }
1353 else {
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);
1357 }
1358 }
1359 }
1360
1361 _moo_median_filter2_model_flat_fibre(rec_fluxes, 10);
1362 {
1363 char *testmodelp3name =
1364 cpl_sprintf("%s_%d_model_step3.txt",
1365 moo_detector_get_extname(det->type, det->ntas),
1366 numfib);
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++) {
1371 int rej;
1372 double v = cpl_image_get(rec_fluxes, x, y, &rej);
1373 fprintf(test5, "%d %d %f\n", x, y - nstep_lo - 1, v);
1374 }
1375 }
1376 cpl_free(testmodelp3name);
1377 fclose(test5);
1378 testmodelp3name =
1379 cpl_sprintf("%s_%d_model_step3.fits",
1380 moo_detector_get_extname(det->type, det->ntas),
1381 numfib);
1382 cpl_image_save(rec_fluxes, testmodelp3name, CPL_TYPE_FLOAT, NULL,
1383 CPL_IO_CREATE);
1384 cpl_free(testmodelp3name);
1385 }
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);
1393
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);
1399
1400 fclose(test4);
1401 cpl_free(testmodel2name);
1402
1403 cpl_free(testname);
1404 fclose(test);
1405
1406 fclose(testp2);
1407 cpl_free(testmodelp2name);
1408 }
1409#endif
1410
1411/* The following addresses an issue with the gcc9 compiler series prior to
1412 * gcc 9.3. These compiler versions require that the variable '__func__'
1413 * appears in the data sharing clause if default(none) is used. However
1414 * adding it to the data sharing clause breaks older compiler versions where
1415 * these variables are pre-determined as shared.
1416 * This was fixed in gcc 9.3 and OpenMP 5.1.
1417 */
1418#ifdef _OPENMP
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)
1423#else
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)
1427#endif
1428 {
1429#pragma omp for
1430#endif
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);
1434 if (h == 0) {
1435 double *data = cpl_image_get_data(rec_fluxes);
1436 for (int j = 0; j < nx * ny; j++) {
1437 data[j] = NAN;
1438 }
1439 }
1440 else {
1441 /*_moo_extract_error_model_flat_fibre(image, i, centroids, nstep_lo, nstep_up,
1442 step, rec_fluxes,nb_fibres);
1443 {
1444 char* testname = cpl_sprintf("error_%d.fits",i);
1445 cpl_image_save(rec_fluxes,testname,CPL_TYPE_DOUBLE,header,CPL_IO_CREATE);
1446 cpl_free(testname);
1447 }*/
1448 /*_moo_rectify_fibre(image, i, centroids, nstep_lo, nstep_up,
1449 step, rec_fluxes);*/
1450 cpl_msg_info(__func__, "Modelling profiles %s: fibre %d",
1451 extname, i);
1452 _moo_model_fibre(image, i, centroids, fwlo, fwup, nstep_lo,
1453 nstep_up, step, rec_fluxes, win_hsize, xstep);
1454
1455 _moo_median_filter2_model_flat_fibre(rec_fluxes, 10);
1456
1457 //_moo_fit_model_flat_fibre(rec_fluxes, deg_poly);
1458 //_moo_smooth_model_flat_fibre(rec_fluxes, 1600);
1459
1460 //_moo_median_filter_model_flat_fibre(rec_fluxes,10,5);
1461 //_moo_savgol_model_flat_fibre(rec_fluxes,3,21);
1462 }
1463 imgtab[i - 1] = rec_fluxes;
1464 }
1465#ifdef _OPENMP
1466 }
1467#endif
1468
1469 cpl_imagelist *cube = cpl_imagelist_new();
1470 for (int i = 0; i < nb_fibres; i++) {
1471 cpl_imagelist_set(cube, imgtab[i], i);
1472 }
1473 cpl_free(imgtab);
1474 res->cube = cube;
1475 moo_psf_single_set_cube_ref(res, nstep_lo + 1, step);
1476#if MOO_DEBUG_MODEL_FLAT
1477 {
1478 double crpix2, cd2_2;
1479 moo_psf_single_get_cube_ref(res, &crpix2, &cd2_2);
1480
1481 char *testname =
1482 cpl_sprintf("%s_neg_cube.txt",
1483 moo_detector_get_extname(det->type, det->ntas));
1484 FILE *test = fopen(testname, "w");
1485 fprintf(test, "#fibnum x y yc flux\n");
1486
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);
1490
1491 for (int y = 1; y <= cy; y++) {
1492 int rej = 0;
1493
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;
1498
1499 double flux = cpl_image_get(img, x, y, &rej);
1500 if (flux < 0) {
1501 fprintf(test, "%d %d %f %d %f\n", f + 1, x, yc, y,
1502 flux);
1503 }
1504 }
1505 }
1506 }
1507 }
1508 fclose(test);
1509 cpl_free(testname);
1510 }
1511#endif
1512 // dump error from the state
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);
1516 // Recover from the error(s) (Reset to prestate))
1517 cpl_errorstate_set(prestate);
1518 }
1519 return res;
1520}
1521
1522/*----------------------------------------------------------------------------*/
1544/*----------------------------------------------------------------------------*/
1545moo_psf *
1546moo_model_flat(moo_det *det,
1547 moo_loc *loc,
1548 moo_model_flat_params *params,
1549 const char *filename)
1550{
1551 moo_psf *result = NULL;
1552
1553 cpl_errorstate prestate = cpl_errorstate_get();
1554
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);
1558
1559 int badpix_level =
1561 // int badpix_level = MOO_BADPIX_OUTSIDE_DATA_RANGE| MOO_BADPIX_COSMETIC;
1562
1563 cpl_msg_info(__func__, "Modelizing flat using oversamp %d",
1564 params->oversamp);
1565 result = moo_psf_create(filename);
1566
1567 cpl_table *fibres_table = moo_loc_get_fibre_table(loc);
1568 cpl_ensure(fibres_table != NULL, CPL_ERROR_ILLEGAL_INPUT, NULL);
1569
1570 cpl_msg_indent_more();
1571 for (int i = 1; i <= 2; i++) {
1572 cpl_table_unselect_all(fibres_table);
1573 int fibname = i;
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);
1577 const int *health =
1578 cpl_table_get_data_int_const(selected, MOO_FIBRES_TABLE_HEALTH);
1579
1580 for (int j = 0; j < 3; j++) {
1581 moo_single *det_single =
1582 moo_det_load_single(det, j, i, badpix_level);
1583 moo_loc_single *loc_single = moo_loc_get_single(loc, j, i);
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);
1590 moo_psf_add_single(result, psf_single, j, i);
1591 }
1592 }
1593 cpl_table_delete(selected);
1594 }
1595 cpl_msg_indent_less();
1596
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);
1600 // Recover from the error(s) (Reset to prestate))
1601 cpl_errorstate_set(prestate);
1602 }
1603 return result;
1604}
#define MOO_BADPIX_OUTSIDE_DATA_RANGE
Definition: moo_badpix.h:62
#define MOO_BADPIX_HOT
Definition: moo_badpix.h:51
#define MOO_BADPIX_COSMETIC
Definition: moo_badpix.h:57
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.
Definition: moo_det.c:197
const char * moo_detector_get_extname(moo_detector_type type, int ntas)
Get the extension name of a detector.
Definition: moo_detector.c:137
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.
Definition: moo_loc.c:155
cpl_table * moo_loc_get_fibre_table(moo_loc *self)
Get the FIBRE TABLE in LOC.
Definition: moo_loc.c:306
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.
Definition: moo_psf.c:228
moo_psf * moo_psf_create(const char *filename)
Create a new empty PSF filename.
Definition: moo_psf.c:82
hdrl_image * moo_single_get_image(moo_single *self)
Get the IMAGE part (DATA,ERR) of single DET.
Definition: moo_single.c:330
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.
Definition: moo_utils.c:1656
double moo_vector_get_min(const cpl_vector *v, int *flags)
Find minimum values in a vector using flags.
Definition: moo_utils.c:979
double moo_vector_get_max(const cpl_vector *v, int *flags)
Find maximum values in a vector using flags.
Definition: moo_utils.c:1003
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.
Definition: moo_utils.c:880
cpl_vector * moo_median_filter(cpl_vector *v, int winhsize)
Apply a median filter to a vector.
Definition: moo_utils.c:1730
double moo_vector_get_percentile(cpl_vector *v, double f)
Get percentile of input vector.
Definition: moo_utils.c:1263