MOONS Pipeline Reference Manual 0.13.1
moo_utils.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
28#include "math.h"
29#include <cpl.h>
30#include <gsl/gsl_multifit.h>
31#include <gsl/gsl_statistics.h>
32#include <gsl/gsl_errno.h>
33#include <gsl/gsl_spline.h>
34#include <gsl/gsl_fit.h>
35#include <gsl/gsl_sort.h>
36
37#include "moo_badpix.h"
38#include "moo_utils.h"
39#include "sc_skycorr.h"
40#include "hdrl_types.h"
41
42#include <string.h>
43
44/*----------------------------------------------------------------------------*/
48/*----------------------------------------------------------------------------*/
49
52static const char *const _moons_license =
53 "This file is part of the MOONS Instrument Pipeline\n"
54 "Copyright (C) 2015-2016 European Southern Observatory\n"
55 "\n"
56 "This program is free software; you can redistribute it and/or modify\n"
57 "it under the terms of the GNU General Public License as published by\n"
58 "the Free Software Foundation; either version 2 of the License, or\n"
59 "(at your option) any later version.\n"
60 "\n"
61 "This program is distributed in the hope that it will be useful,\n"
62 "but WITHOUT ANY WARRANTY; without even the implied warranty of\n"
63 "MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the\n"
64 "GNU General Public License for more details.\n"
65 "\n"
66 "You should have received a copy of the GNU General Public License\n"
67 "along with this program; if not, write to the Free Software\n"
68 "Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,\n"
69 "MA 02110-1301 USA.";
70
71/*----------------------------------------------------------------------------*/
79/*----------------------------------------------------------------------------*/
80const char *
82{
83 return _moons_license;
84}
85
86/*----------------------------------------------------------------------------*/
98/*----------------------------------------------------------------------------*/
99cpl_image *
100moo_compute_sigma_map(cpl_imagelist *list, cpl_imagelist *qlist, cpl_image *img)
101{
102 cpl_ensure(list != NULL, CPL_ERROR_NULL_INPUT, NULL);
103 cpl_ensure(qlist != NULL, CPL_ERROR_NULL_INPUT, NULL);
104 cpl_ensure(img != NULL, CPL_ERROR_NULL_INPUT, NULL);
105
106 int size = cpl_imagelist_get_size(list);
107 int qsize = cpl_imagelist_get_size(qlist);
108 cpl_ensure(size > 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
109 cpl_ensure(size == qsize, CPL_ERROR_ILLEGAL_INPUT, NULL);
110
111 int nx = cpl_image_get_size_x(img);
112 int ny = cpl_image_get_size_y(img);
113 const double *med_data = cpl_image_get_data_double_const(img);
114 int i, j;
115 cpl_image *res = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
116 double *res_data = cpl_image_get_data_double(res);
117
118 for (i = 0; i < nx * ny; i++) {
119 double med = med_data[i];
120
121 double val = 0;
122 int nbcomb = 0;
123 for (j = 0; j < size; j++) {
124 const cpl_image *im = cpl_imagelist_get_const(list, j);
125 const cpl_image *qim = cpl_imagelist_get_const(qlist, j);
126 const double *im_data = cpl_image_get_data_double_const(im);
127 const int *qim_data = cpl_image_get_data_int_const(qim);
128
129 if (qim_data[i] == MOO_BADPIX_GOOD) {
130 double sig = im_data[i] - med;
131 val += sig * sig;
132 nbcomb++;
133 }
134 }
135 if (nbcomb > 0) {
136 res_data[i] = sqrt(val / (size - 1));
137 }
138 else {
139 res_data[i] = NAN;
140 }
141 }
142 cpl_image_set_bpm(res, cpl_mask_duplicate(cpl_image_get_bpm(img)));
143
144 return res;
145}
146
147/*----------------------------------------------------------------------------*/
161/*----------------------------------------------------------------------------*/
162cpl_mask *
163moo_kappa_sigma_clipping(cpl_image *sigma_img,
164 int niter,
165 double kappa,
166 double cdiff,
167 double maxfrac)
168{
169 cpl_mask *res = NULL;
170 cpl_ensure(sigma_img != NULL, CPL_ERROR_NULL_INPUT, NULL);
171 cpl_ensure(niter >= 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
172 cpl_ensure(kappa >= 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
173
174 int iter, i;
175
176 int nx = cpl_image_get_size_x(sigma_img);
177 int ny = cpl_image_get_size_y(sigma_img);
178 cpl_mask *mask = cpl_image_get_bpm(sigma_img);
179
180 res = cpl_mask_new(nx, ny);
181
182 cpl_binary *orig_mask_data = cpl_mask_get_data(mask);
183 cpl_binary *res_mask_data = cpl_mask_get_data(res);
184
185 double *sigma_data = cpl_image_get_data_double(sigma_img);
186
187 int initbad = cpl_mask_count(mask);
188 int initgood = nx * ny - initbad;
189 double maxrejected = maxfrac * initgood;
190
191 double sig = 0;
192 double median = cpl_image_get_median_dev(sigma_img, &sig);
193
194 for (iter = 0; iter < niter; iter++) {
195 cpl_msg_info(__func__,
196 "Iteration %d: median value = %f stdev value = %f",
197 iter + 1, median, sig);
198 int nbrej = 0;
199
200 double ksig = kappa * sig;
201 for (i = 0; i < nx * ny; i++) {
202 if (orig_mask_data[i] == CPL_BINARY_0 &&
203 fabs(sigma_data[i] - median) > ksig) {
204 nbrej++;
205 if (nbrej < maxrejected) {
206 res_mask_data[i] = CPL_BINARY_1;
207 orig_mask_data[i] = CPL_BINARY_1;
208 }
209 else {
210 cpl_msg_info(
211 __func__,
212 "Number of rejected pixel reached maximum value %d/%d",
213 nbrej, initgood);
214 break;
215 }
216 }
217 }
218 if (nbrej == 0) {
219 break;
220 }
221 cpl_msg_info(__func__, "Number of detected bad pixels %d", nbrej);
222 double sig2 = 0;
223 median = cpl_image_get_median_dev(sigma_img, &sig2);
224
225 if ((sig - sig2) > cdiff) {
226 sig = sig2;
227 }
228 else {
229 cpl_msg_info(__func__,
230 "Difference in sigma = %f reached minimum value = %f",
231 (sig - sig2), cdiff);
232 break;
233 }
234 }
235
236 return res;
237}
238
239
240/*----------------------------------------------------------------------------*/
252/*----------------------------------------------------------------------------*/
253cpl_error_code
254moo_barycenter_fit(cpl_bivector *points, double *center, double *width)
255{
256 cpl_ensure_code(points != NULL, CPL_ERROR_NULL_INPUT);
257 cpl_ensure_code(center != NULL, CPL_ERROR_NULL_INPUT);
258 cpl_ensure_code(width != NULL, CPL_ERROR_NULL_INPUT);
259
260 int size = cpl_bivector_get_size(points);
261
262 cpl_ensure_code(size > 0, CPL_ERROR_ILLEGAL_INPUT);
263
264 cpl_vector *x = cpl_bivector_get_x(points);
265 cpl_vector *y = cpl_bivector_get_y(points);
266
267 double x1 = cpl_vector_get(x, 0);
268 double x2 = cpl_vector_get(x, size - 1);
269
270 double res = 0.0;
271 double sum = 0.0;
272 int i;
273 for (i = 0; i < size; i++) {
274 double p = cpl_vector_get(x, i);
275 double f = cpl_vector_get(y, i);
276 if (f > 0) {
277 res += p * f;
278 sum += f;
279 }
280 }
281 if (sum == 0) {
282 return CPL_ERROR_DIVISION_BY_ZERO;
283 }
284 res /= sum;
285 *center = res;
286 *width = x2 - x1;
287
288 return CPL_ERROR_NONE;
289}
290
291/*----------------------------------------------------------------------------*/
307/*----------------------------------------------------------------------------*/
308cpl_error_code
309moo_gaussian_fit(cpl_bivector *points,
310 cpl_fit_mode fit_pars,
311 double *center,
312 double *width,
313 double *background,
314 double *area)
315{
316 cpl_ensure_code(points != NULL, CPL_ERROR_NULL_INPUT);
317 cpl_ensure_code(center != NULL, CPL_ERROR_NULL_INPUT);
318 cpl_ensure_code(width != NULL, CPL_ERROR_NULL_INPUT);
319 cpl_ensure_code(background != NULL, CPL_ERROR_NULL_INPUT);
320 cpl_ensure_code(area != NULL, CPL_ERROR_NULL_INPUT);
321
322 int size = cpl_bivector_get_size(points);
323
324 cpl_ensure_code(size > 0, CPL_ERROR_ILLEGAL_INPUT);
325
326 cpl_vector *x = cpl_bivector_get_x(points);
327 cpl_vector *y = cpl_bivector_get_y(points);
328
329 cpl_error_code status =
330 cpl_vector_fit_gaussian(x, NULL, y, NULL, fit_pars, center, width, area,
331 background, NULL, NULL, NULL);
332
333 return status;
334}
335
336/*----------------------------------------------------------------------------*/
346/*----------------------------------------------------------------------------*/
347double
348moo_gaussian_eval(double x, double x0, double sigma, double offset, double area)
349{
350 double res;
351
352 res = area / sqrt(CPL_MATH_2PI * sigma * sigma) *
353 exp(-(x - x0) * (x - x0) / (2 * sigma * sigma)) +
354 offset;
355
356 return res;
357}
358
359/*----------------------------------------------------------------------------*/
371/*----------------------------------------------------------------------------*/
372cpl_error_code
374 double x0,
375 double sigma,
376 double offset,
377 double area,
378 double *x1,
379 double *x2)
380{
381 double n = log((y - offset) * (sqrt(CPL_MATH_2PI * sigma * sigma) / area)) *
382 (-2 * sigma * sigma);
383
384 *x1 = x0 - sqrt(n);
385 *x2 = x0 + sqrt(n);
386
387 return CPL_ERROR_NONE;
388}
389
390static double
391_moo_interpolate_x(cpl_bivector *points, int i, double val)
392{
393 cpl_ensure_code(points != NULL, CPL_ERROR_NULL_INPUT);
394 cpl_ensure_code(i >= 0, CPL_ERROR_NULL_INPUT);
395
396 int size = cpl_bivector_get_size(points);
397
398
399 cpl_ensure_code(size > i + 1, CPL_ERROR_NULL_INPUT);
400
401 double *y = cpl_bivector_get_y_data(points);
402 double *x = cpl_bivector_get_x_data(points);
403
404 double x1 = x[i];
405 double y1 = y[i];
406 double x2 = x[i + 1];
407 double y2 = y[i + 1];
408
409 return x1 + (val - y1) / (y2 - y1) * (x2 - x1);
410}
411
412/*----------------------------------------------------------------------------*/
421/*----------------------------------------------------------------------------*/
422cpl_error_code
423moo_find_threshold_limits(cpl_bivector *points,
424 double threshold,
425 double *xmin,
426 double *xmax)
427{
428 cpl_ensure_code(points != NULL, CPL_ERROR_NULL_INPUT);
429 cpl_ensure_code(xmin != NULL, CPL_ERROR_NULL_INPUT);
430 cpl_ensure_code(xmax != NULL, CPL_ERROR_NULL_INPUT);
431
432 int size = cpl_bivector_get_size(points);
433
434 cpl_ensure_code(size > 1, CPL_ERROR_NULL_INPUT);
435
436 double *y = cpl_bivector_get_y_data(points);
437 int i;
438
439 for (i = 1; i < size; i++) {
440 if (y[i] > threshold) {
441 *xmin = _moo_interpolate_x(points, i - 1, threshold);
442 break;
443 }
444 }
445
446 for (i = size - 2; i >= 0; i++) {
447 if (y[i] > threshold) {
448 *xmax = _moo_interpolate_x(points, i, threshold);
449 break;
450 }
451 }
452 return CPL_ERROR_NONE;
453}
454
455
456/*----------------------------------------------------------------------------*/
464/*----------------------------------------------------------------------------*/
465static cpl_error_code
466_moo_tchebychev_transform(cpl_vector *posv, double min, double max)
467{
468 int i;
469 double a, b;
470
471 cpl_ensure_code(posv != NULL, CPL_ERROR_NULL_INPUT);
472 int size = cpl_vector_get_size(posv);
473 cpl_ensure_code(size > 0, CPL_ERROR_ILLEGAL_INPUT);
474 cpl_ensure_code(min < max, CPL_ERROR_ILLEGAL_INPUT);
475
476 double *pos = cpl_vector_get_data(posv);
477
478 a = 2 / (max - min);
479 b = 1 - 2 * max / (max - min);
480
481 for (i = 0; i < size; i++) {
482 double res = pos[i] * a + b;
483 if (res < -1.) {
484 pos[i] = -1.;
485 }
486 else if (res > 1.) {
487 pos[i] = 1.;
488 }
489 else {
490 pos[i] = res;
491 }
492 }
493 return CPL_ERROR_NONE;
494}
495
496static cpl_error_code
497_moo_tchebychev_reverse_transform(cpl_vector *posv, double min, double max)
498{
499 double a, b;
500 int i;
501
502 cpl_ensure_code(min < max, CPL_ERROR_ILLEGAL_INPUT);
503
504 int size = cpl_vector_get_size(posv);
505
506 a = 2 / (max - min);
507 b = 1 - 2 * max / (max - min);
508
509 double *pos = cpl_vector_get_data(posv);
510
511 for (i = 0; i < size; i++) {
512 double res = (pos[i] - b) / a;
513 pos[i] = res;
514 }
515
516 return CPL_ERROR_NONE;
517}
518/*----------------------------------------------------------------------------*/
529/*----------------------------------------------------------------------------*/
530static cpl_vector *
531_moo_tchebychev_poly_eval(int n, double X)
532{
533 cpl_vector *result = NULL;
534
535 cpl_ensure(n >= 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
536
537 result = cpl_vector_new(n + 1);
538
539 cpl_vector_set(result, 0, 1.0);
540
541 if (n > 0) {
542 int indice = 2;
543 cpl_vector_set(result, 1, X);
544 while (indice <= n) {
545 double T_indice = 0.0;
546 double T_indice_1 = 0.0;
547 double T_indice_2 = 0.0;
548
549 T_indice_1 = cpl_vector_get(result, indice - 1);
550 T_indice_2 = cpl_vector_get(result, indice - 2);
551 T_indice = 2 * X * T_indice_1 - T_indice_2;
552 cpl_vector_set(result, indice, T_indice);
553 indice++;
554 }
555 }
556 return result;
557}
558
559/*----------------------------------------------------------------------------*/
571/*----------------------------------------------------------------------------*/
572moo_tcheby_polynomial *
573moo_tcheby_polynomial_fit(cpl_bivector *data,
574 int degree,
575 double xmin,
576 double xmax)
577{
578 cpl_ensure(data != NULL, CPL_ERROR_NULL_INPUT, NULL);
579 cpl_ensure(degree > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
580
581 int i, j;
582 double chisq;
583 gsl_matrix *X = NULL, *cov = NULL;
584 gsl_vector *y = NULL, *c = NULL;
585 gsl_multifit_linear_workspace *work = NULL;
586 cpl_size pows[1];
587
588 cpl_ensure(data != NULL, CPL_ERROR_NULL_INPUT, NULL);
589
590 int size = cpl_bivector_get_size(data);
591
592 cpl_ensure(size > degree, CPL_ERROR_ILLEGAL_INPUT, NULL);
593
594 moo_tcheby_polynomial *res = cpl_calloc(1, sizeof(moo_tcheby_polynomial));
595
596 res->solution = cpl_polynomial_new(1);
597
598 cpl_vector *vx = cpl_bivector_get_x(data);
599 cpl_vector *vy = cpl_bivector_get_y(data);
600
601 res->xmax = xmax;
602 res->xmin = xmin;
603 res->ymin = cpl_vector_get_min(vy);
604 res->ymax = cpl_vector_get_max(vy);
605
606 res->degree = degree;
607
608 _moo_tchebychev_transform(vx, res->xmin, res->xmax);
609 _moo_tchebychev_transform(vy, res->ymin, res->ymax);
610 double *xdata = cpl_bivector_get_x_data(data);
611 double *ydata = cpl_bivector_get_y_data(data);
612
613 int nbcoeffs = degree + 1;
614
615 X = gsl_matrix_alloc(size, nbcoeffs);
616 y = gsl_vector_alloc(size);
617 c = gsl_vector_alloc(nbcoeffs);
618 cov = gsl_matrix_alloc(nbcoeffs, nbcoeffs);
619
620 for (i = 0; i < size; i++) {
621 for (j = 0; j < degree + 1; j++) {
622 double Tj = cos(j * acos(xdata[i]));
623 gsl_matrix_set(X, i, j, Tj);
624 }
625 gsl_vector_set(y, i, ydata[i]);
626 }
627
628 work = gsl_multifit_linear_alloc(size, nbcoeffs);
629 gsl_multifit_linear(X, y, c, cov, &chisq, work);
630
631 for (j = 0; j < degree + 1; j++) {
632 pows[0] = j;
633 cpl_polynomial_set_coeff(res->solution, pows, gsl_vector_get(c, j));
634 }
635
636 _moo_tchebychev_reverse_transform(vy, res->ymin, res->ymax);
637 _moo_tchebychev_reverse_transform(vx, res->xmin, res->xmax);
638 gsl_multifit_linear_free(work);
639 gsl_matrix_free(X);
640 gsl_vector_free(y);
641 gsl_vector_free(c);
642 gsl_matrix_free(cov);
643
644 return res;
645}
646
647/*----------------------------------------------------------------------------*/
665/*----------------------------------------------------------------------------*/
666moo_tcheby2d_polynomial *
668 int xdegree,
669 double xmin,
670 double xmax,
671 cpl_vector *in_y,
672 int ydegree,
673 double ymin,
674 double ymax,
675 cpl_vector *in_l,
676 double lmin,
677 double lmax)
678{
679 cpl_ensure(in_x != NULL, CPL_ERROR_NULL_INPUT, NULL);
680 cpl_ensure(xdegree > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
681 cpl_ensure(in_y != NULL, CPL_ERROR_NULL_INPUT, NULL);
682 cpl_ensure(ydegree > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
683 cpl_ensure(in_l != NULL, CPL_ERROR_NULL_INPUT, NULL);
684
685 int i, j, k;
686 double chisq;
687 gsl_matrix *X = NULL, *cov = NULL;
688 gsl_vector *y = NULL, *c = NULL;
689 gsl_multifit_linear_workspace *work = NULL;
690 cpl_size pows[2];
691
692 moo_tcheby2d_polynomial *res =
693 cpl_calloc(1, sizeof(moo_tcheby2d_polynomial));
694
695 res->solution = cpl_polynomial_new(2);
696
697 int size = cpl_vector_get_size(in_x);
698
699 res->xmin = xmin;
700 res->xmax = xmax;
701 res->ymin = ymin;
702 res->ymax = ymax;
703 res->lmin = lmin;
704 res->lmax = lmax;
705 res->degree_x = xdegree;
706 res->degree_y = ydegree;
707
708 cpl_vector *vx = cpl_vector_duplicate(in_x);
709 cpl_vector *vy = cpl_vector_duplicate(in_y);
710 cpl_vector *vl = cpl_vector_duplicate(in_l);
711
712 _moo_tchebychev_transform(vx, res->xmin, res->xmax);
713 _moo_tchebychev_transform(vy, res->ymin, res->ymax);
714 _moo_tchebychev_transform(vl, res->lmin, res->lmax);
715
716 double *xdata = cpl_vector_get_data(vx);
717 double *ydata = cpl_vector_get_data(vy);
718 double *ldata = cpl_vector_get_data(vl);
719
720 int nbcoeffs = (xdegree + 1) * (ydegree + 1);
721
722 X = gsl_matrix_alloc(size, nbcoeffs);
723 y = gsl_vector_alloc(size);
724 c = gsl_vector_alloc(nbcoeffs);
725 cov = gsl_matrix_alloc(nbcoeffs, nbcoeffs);
726
727 for (i = 0; i < size; i++) {
728 for (j = 0; j < xdegree + 1; j++) {
729 for (k = 0; k < ydegree + 1; k++) {
730 double px = cos(j * acos(xdata[i]));
731 double py = cos(k * acos(ydata[i]));
732 gsl_matrix_set(X, i, k + j * (ydegree + 1), px * py);
733 }
734 }
735 gsl_vector_set(y, i, ldata[i]);
736 }
737 work = gsl_multifit_linear_alloc(size, nbcoeffs);
738 gsl_multifit_linear(X, y, c, cov, &chisq, work);
739
740 for (j = 0; j < xdegree + 1; j++) {
741 for (k = 0; k < ydegree + 1; k++) {
742 pows[0] = k;
743 pows[1] = j;
744 cpl_polynomial_set_coeff(res->solution, pows,
745 gsl_vector_get(c, k + j * (ydegree + 1)));
746 }
747 }
748
749 gsl_multifit_linear_free(work);
750 gsl_matrix_free(X);
751 gsl_vector_free(y);
752 gsl_vector_free(c);
753 gsl_matrix_free(cov);
754
755 cpl_vector_delete(vx);
756 cpl_vector_delete(vy);
757 cpl_vector_delete(vl);
758
759 return res;
760}
761
762void
763moo_tcheby_polynomial_delete(moo_tcheby_polynomial *self)
764{
765 if (self != NULL) {
766 if (self->solution != NULL) {
767 cpl_polynomial_delete(self->solution);
768 }
769 cpl_free(self);
770 }
771}
772
773void
774moo_tcheby2d_polynomial_delete(moo_tcheby2d_polynomial *self)
775{
776 if (self != NULL) {
777 if (self->solution != NULL) {
778 cpl_polynomial_delete(self->solution);
779 }
780 cpl_free(self);
781 }
782}
783
784double
785moo_tcheby_polynomial_eval(moo_tcheby_polynomial *self, double x)
786{
787 cpl_size pows[1];
788 double yval = 0;
789
790 cpl_vector *vx = cpl_vector_new(1);
791 cpl_vector *vy = cpl_vector_new(1);
792
793 cpl_vector_set(vx, 0, x);
794
795 _moo_tchebychev_transform(vx, self->xmin, self->xmax);
796 double tx = cpl_vector_get(vx, 0);
797
798 cpl_vector *tcheb = _moo_tchebychev_poly_eval(self->degree, tx);
799 for (int j = 0; j < self->degree + 1; j++) {
800 pows[0] = j;
801 double tval = cpl_vector_get(tcheb, j);
802 double coef = cpl_polynomial_get_coeff(self->solution, pows);
803 yval += coef * tval;
804 }
805 cpl_vector_set(vy, 0, yval);
806
807 _moo_tchebychev_reverse_transform(vy, self->ymin, self->ymax);
808
809 yval = cpl_vector_get(vy, 0);
810
811 cpl_vector_delete(tcheb);
812 cpl_vector_delete(vx);
813 cpl_vector_delete(vy);
814
815 return yval;
816}
817
818double
819moo_tcheby2d_polynomial_eval(moo_tcheby2d_polynomial *self, double x, double y)
820{
821 cpl_size pows[2];
822 double yval = 0;
823
824 cpl_vector *vx = cpl_vector_new(1);
825 cpl_vector *vy = cpl_vector_new(1);
826 cpl_vector *vl = cpl_vector_new(1);
827
828 cpl_vector_set(vx, 0, x);
829 cpl_vector_set(vy, 0, y);
830
831 _moo_tchebychev_transform(vx, self->xmin, self->xmax);
832 _moo_tchebychev_transform(vy, self->ymin, self->ymax);
833 double tx = cpl_vector_get(vx, 0);
834 double ty = cpl_vector_get(vy, 0);
835
836 cpl_vector *tchebx = _moo_tchebychev_poly_eval(self->degree_x, tx);
837 cpl_vector *tcheby = _moo_tchebychev_poly_eval(self->degree_y, ty);
838
839 for (int j = 0; j < self->degree_x + 1; j++) {
840 for (int k = 0; k < self->degree_y + 1; k++) {
841 pows[0] = k;
842 pows[1] = j;
843 double tvalx = cpl_vector_get(tchebx, j);
844 double tvaly = cpl_vector_get(tcheby, k);
845 double coef = cpl_polynomial_get_coeff(self->solution, pows);
846 yval += coef * tvalx * tvaly;
847 }
848 }
849
850 cpl_vector_set(vl, 0, yval);
851
852 _moo_tchebychev_reverse_transform(vl, self->lmin, self->lmax);
853
854 yval = cpl_vector_get(vl, 0);
855
856 cpl_vector_delete(tchebx);
857 cpl_vector_delete(tcheby);
858 cpl_vector_delete(vx);
859 cpl_vector_delete(vy);
860 cpl_vector_delete(vl);
861 return yval;
862}
863/*----------------------------------------------------------------------------*/
878/*----------------------------------------------------------------------------*/
879cpl_error_code
880moo_tchebychev_fit(cpl_bivector *data,
881 int *flag,
882 int degree,
883 double xmin,
884 double xmax,
885 double ymin,
886 double ymax)
887{
888 cpl_ensure_code(data != NULL, CPL_ERROR_NULL_INPUT);
889 cpl_ensure_code(flag != NULL, CPL_ERROR_NULL_INPUT);
890 cpl_ensure_code(degree > 0, CPL_ERROR_ILLEGAL_INPUT);
891 int i, j;
892 double chisq;
893 gsl_matrix *X = NULL, *cov = NULL;
894 gsl_vector *y = NULL, *c = NULL;
895 gsl_multifit_linear_workspace *work = NULL;
896 cpl_size pows[1];
897 cpl_polynomial *result = cpl_polynomial_new(1);
898
899 cpl_vector *vx = cpl_bivector_get_x(data);
900 cpl_vector *vy = cpl_bivector_get_y(data);
901
902 _moo_tchebychev_transform(vx, xmin, xmax);
903 _moo_tchebychev_transform(vy, ymin, ymax);
904
905
906 double *xdata = cpl_bivector_get_x_data(data);
907 double *ydata = cpl_bivector_get_y_data(data);
908
909 int nbcoeffs = degree + 1;
910
911 int nb_valid = 0;
912 for (i = xmin - 1; i < xmax; i++) {
913 if (flag[i] == 0) {
914 nb_valid++;
915 }
916 }
917
918 X = gsl_matrix_alloc(nb_valid, nbcoeffs);
919 y = gsl_vector_alloc(nb_valid);
920 c = gsl_vector_alloc(nbcoeffs);
921 cov = gsl_matrix_alloc(nbcoeffs, nbcoeffs);
922
923 int id = 0;
924 for (i = xmin - 1; i < xmax; i++) {
925 if (flag[i] == 0) {
926 for (j = 0; j < degree + 1; j++) {
927 double Tj = cos(j * acos(xdata[i]));
928 gsl_matrix_set(X, id, j, Tj);
929 }
930 gsl_vector_set(y, id, ydata[i]);
931 id++;
932 }
933 }
934
935 work = gsl_multifit_linear_alloc(nb_valid, nbcoeffs);
936 gsl_multifit_linear(X, y, c, cov, &chisq, work);
937
938 for (j = 0; j < degree + 1; j++) {
939 pows[0] = j;
940 cpl_polynomial_set_coeff(result, pows, gsl_vector_get(c, j));
941 }
942
943 for (i = xmin - 1; i < xmax; i++) {
944 double yval = 0;
945 double xval = xdata[i];
946 cpl_vector *tcheb = _moo_tchebychev_poly_eval(degree, xval);
947 for (j = 0; j < degree + 1; j++) {
948 pows[0] = j;
949 double tval = cpl_vector_get(tcheb, j);
950 double coef = cpl_polynomial_get_coeff(result, pows);
951 yval += coef * tval;
952 ydata[i] = yval;
953 }
954 cpl_vector_delete(tcheb);
955 }
956
957 _moo_tchebychev_reverse_transform(vx, xmin, xmax);
958 _moo_tchebychev_reverse_transform(vy, ymin, ymax);
959
960 gsl_multifit_linear_free(work);
961 gsl_matrix_free(X);
962 gsl_vector_free(y);
963 gsl_vector_free(c);
964 gsl_matrix_free(cov);
965 cpl_polynomial_delete(result);
966
967 return CPL_ERROR_NONE;
968}
969
970/*----------------------------------------------------------------------------*/
977/*----------------------------------------------------------------------------*/
978double
979moo_vector_get_min(const cpl_vector *v, int *flags)
980{
981 int size = cpl_vector_get_size(v);
982
983 double res = DBL_MAX;
984
985 for (int i = 0; i < size; i++) {
986 if (flags[i] == 0) {
987 double val = cpl_vector_get(v, i);
988 res = CPL_MIN(res, val);
989 }
990 }
991 return res;
992}
993
994/*----------------------------------------------------------------------------*/
1001/*----------------------------------------------------------------------------*/
1002double
1003moo_vector_get_max(const cpl_vector *v, int *flags)
1004{
1005 int size = cpl_vector_get_size(v);
1006
1007 double res = DBL_MIN;
1008
1009
1010 for (int i = 0; i < size; i++) {
1011 if (flags[i] == 0) {
1012 double val = cpl_vector_get(v, i);
1013 res = CPL_MAX(res, val);
1014 }
1015 }
1016 return res;
1017}
1018
1019cpl_vector *
1020moo_hpfilter(cpl_vector *Yv, double S)
1021{
1022 int N = cpl_vector_get_size(Yv);
1023 double *Y = cpl_vector_get_data(Yv);
1024
1025 cpl_vector *Tv = cpl_vector_new(N);
1026 double *T = cpl_vector_get_data(Tv);
1027 double **V = cpl_calloc(N, sizeof(double *));
1028
1029 for (int i = 0; i < N; i++) {
1030 V[i] = cpl_calloc(3, sizeof(double));
1031 }
1032 double *D = cpl_calloc(N, sizeof(double));
1033 double SS = 0.0, M1 = 0.0, M2 = 0.0, V11 = 0.0, V12 = 0.0, V22 = 0.0,
1034 X = 0.0, Z = 0.0, B11 = 0.0, B12 = 0.0, B22 = 0.0, DET = 0.0,
1035 E1 = 0.0, E2 = 0.0;
1036 int NN = 0, I1 = 0, IB = 0;
1037
1038 if (NN != N || S != SS) {
1039 SS = S;
1040 NN = N;
1041 V11 = 1;
1042 V22 = 1;
1043 V12 = 0;
1044 for (int i = 2; i <= N - 1; i++) {
1045 X = V11;
1046 Z = V12;
1047 V11 = 1.0 / S + 4.0 * (X - Z) + V22;
1048 V12 = 2.0 * X - Z;
1049 V22 = X;
1050 DET = V11 * V22 - V12 * V12;
1051 V[i][0] = V22 / DET;
1052 V[i][1] = -V12 / DET;
1053 V[i][2] = V11 / DET;
1054 X = V11 + 1;
1055 Z = V11;
1056 V11 = V11 - V11 * V11 / X;
1057 V22 = V22 - V12 * V12 / X;
1058 V12 = V12 - Z * V12 / X;
1059 }
1060 }
1061
1062 // this is the forward pass
1063 M1 = Y[1];
1064 M2 = Y[0];
1065
1066 for (int i = 2; i <= N - 1; i++) {
1067 X = M1;
1068 M1 = 2.0 * M1 - M2;
1069 M2 = X;
1070 T[i - 1] = V[i][0] * M1 + V[i][1] * M2;
1071 D[i - 1] = V[i][1] * M1 + V[i][2] * M2;
1072 DET = V[i][0] * V[i][2] - V[i][1] * V[i][1];
1073 V11 = V[i][2] / DET;
1074 V12 = -V[i][1] / DET;
1075 Z = (Y[i] - M1) / (V11 + 1);
1076 M1 = M1 + V11 * Z;
1077 M2 = M2 + V12 * Z;
1078 }
1079 T[N - 1] = M1;
1080 T[N - 2] = M2;
1081
1082 //this is the backward pass
1083 M1 = Y[N - 2];
1084 M2 = Y[N - 1];
1085
1086 for (int i = N - 3; i >= 0; i--) {
1087 I1 = i + 1;
1088 IB = (N - i) - 1;
1089 X = M1;
1090 M1 = 2 * M1 - M2;
1091 M2 = X;
1092 // combine info for y(.lt.i) with info for y(.ge.i)
1093 if (i > 2) {
1094 E1 = V[IB][2] * M2 + V[IB][1] * M1 + T[i];
1095 E2 = V[IB][1] * M2 + V[IB][0] * M1 + D[i];
1096 B11 = V[IB][2] + V[I1][0];
1097 B12 = V[IB][1] + V[I1][1];
1098 B22 = V[IB][0] + V[I1][2];
1099 DET = B11 * B22 - B12 * B12;
1100 T[i] = (-B12 * E1 + B11 * E2) / DET;
1101 }
1102 // end of combining
1103 DET = V[IB][0] * V[IB][2] - V[IB][1] * V[IB][1];
1104 V11 = V[IB][2] / DET;
1105 V12 = -V[IB][1] / DET;
1106 Z = (Y[i] - M1) / (V11 + 1.);
1107 M1 = M1 + V11 * Z;
1108 M2 = M2 + V12 * Z;
1109 }
1110
1111 T[0] = M1;
1112 T[1] = M2;
1113
1114 for (int i = 0; i < N; i++) {
1115 cpl_free(V[i]);
1116 }
1117 cpl_free(V);
1118 cpl_free(D);
1119
1120 return Tv;
1121}
1122/*----------------------------------------------------------------------------*/
1157/*----------------------------------------------------------------------------*/
1158
1159cpl_error_code
1160moo_interpolate_linear(cpl_bivector *fout,
1161 cpl_vector *fout_errs,
1162 int *fout_qual,
1163 const cpl_bivector *fref,
1164 const cpl_vector *fref_errs,
1165 const int *fref_qual)
1166{
1167 const cpl_size m = cpl_bivector_get_size(fref);
1168 const cpl_size n = cpl_bivector_get_size(fout);
1169 const double *xref = cpl_bivector_get_x_data_const(fref);
1170 const double *yref = cpl_bivector_get_y_data_const(fref);
1171 const double *ref_errs = cpl_vector_get_data_const(fref_errs);
1172 double *xout = cpl_bivector_get_x_data(fout);
1173 double *yout = cpl_bivector_get_y_data(fout);
1174 double *out_errs = cpl_vector_get_data(fout_errs);
1175
1176 int qual = 0;
1177 double y0_err = DBL_MAX;
1178 double grad = DBL_MAX; /* Avoid (false) uninit warning */
1179 double grad_err = DBL_MAX; /* Avoid (false) uninit warning */
1180 double y_0 = DBL_MAX; /* Avoid (false) uninit warning */
1181 cpl_size ibelow, iabove;
1182 cpl_size i;
1183
1184
1185 cpl_ensure_code(xref != NULL, cpl_error_get_code());
1186 cpl_ensure_code(yref != NULL, cpl_error_get_code());
1187 cpl_ensure_code(ref_errs != NULL, cpl_error_get_code());
1188 cpl_ensure_code(xout != NULL, cpl_error_get_code());
1189 cpl_ensure_code(yout != NULL, cpl_error_get_code());
1190 cpl_ensure_code(out_errs != NULL, cpl_error_get_code());
1191
1192 /* Upper extrapolation not allowed */
1193 cpl_ensure_code(xout[n - 1] <= xref[m - 1], CPL_ERROR_DATA_NOT_FOUND);
1194
1195 /* Start interpolation from below */
1196 ibelow = cpl_vector_find(cpl_bivector_get_x_const(fref), xout[0]);
1197
1198 if (xout[0] < xref[ibelow]) {
1199 /* Lower extrapolation also not allowed */
1200 cpl_ensure_code(ibelow > 0, CPL_ERROR_DATA_NOT_FOUND);
1201 ibelow--;
1202 }
1203
1204 iabove = ibelow; /* Ensure grad initialization, for 1st interpolation */
1205
1206 /* assert( xref[iabove] <= xout[0] ); */
1207 for (i = 0; i < n; i++) {
1208 /* When possible reuse reference function abscissa points */
1209 if (xref[iabove] < xout[i]) {
1210 /* No, need new points */
1211 while (xref[++iabove] < xout[i])
1212 ;
1213 ibelow = iabove - 1;
1214
1215 /* Verify that the pair of reference abscissa points are valid */
1216 if (xref[iabove] <= xref[ibelow])
1217 break;
1218
1219 qual = fref_qual[ibelow] | fref_qual[iabove];
1220
1221 grad =
1222 (yref[iabove] - yref[ibelow]) / (xref[iabove] - xref[ibelow]);
1223 y_0 = yref[ibelow] - grad * xref[ibelow];
1224
1225 grad_err = (ref_errs[iabove] * ref_errs[iabove] -
1226 ref_errs[ibelow] * ref_errs[ibelow]) /
1227 (xref[iabove] - xref[ibelow]);
1228 y0_err =
1229 ref_errs[ibelow] * ref_errs[ibelow] - grad_err * xref[ibelow];
1230 }
1231
1232 if (xref[iabove] > xout[i]) {
1233 if (xref[ibelow] < xout[i]) {
1234 fout_qual[i] = qual;
1235 out_errs[i] = sqrt(y0_err + grad_err * xout[i]);
1236 yout[i] = y_0 + grad * xout[i];
1237 }
1238 else {
1239 /* assert( xout[i] == xref[ibelow] ); */
1240 yout[i] = yref[ibelow];
1241 out_errs[i] = ref_errs[ibelow];
1242 fout_qual[i] = fref_qual[ibelow];
1243 }
1244 }
1245 else {
1246 yout[i] = yref[iabove];
1247 out_errs[i] = ref_errs[iabove];
1248 fout_qual[i] = fref_qual[iabove];
1249 }
1250 }
1251
1252 return i == n ? CPL_ERROR_NONE : CPL_ERROR_ILLEGAL_INPUT;
1253}
1254
1255/*----------------------------------------------------------------------------*/
1262double
1263moo_vector_get_percentile(cpl_vector *v, double f)
1264{
1265 double res = NAN;
1266
1267 cpl_ensure(v != NULL, CPL_ERROR_NULL_INPUT, res);
1268 cpl_ensure(f >= 0 && f <= 1, CPL_ERROR_ILLEGAL_INPUT, res);
1269
1270 int size = cpl_vector_get_size(v);
1271 int nbnan = 0;
1272
1273 for (int i = 0; i < size; i++) {
1274 double val = cpl_vector_get(v, i);
1275 if (isnan(val)) {
1276 nbnan++;
1277 }
1278 }
1279 cpl_vector *temp = cpl_vector_new(size - nbnan);
1280 int idx = 0;
1281 for (int i = 0; i < size; i++) {
1282 double val = cpl_vector_get(v, i);
1283 if (!isnan(val)) {
1284 cpl_vector_set(temp, idx, val);
1285 idx++;
1286 }
1287 }
1288
1289 cpl_vector_sort(temp, CPL_SORT_ASCENDING);
1290
1291 double *data = cpl_vector_get_data(temp);
1292 int dsize = cpl_vector_get_size(temp);
1293
1294 res = gsl_stats_quantile_from_sorted_data(data, 1, dsize, f);
1295
1296 cpl_vector_delete(temp);
1297
1298 return res;
1299}
1300
1301/*----------------------------------------------------------------------------*/
1308cpl_vector *
1310{
1311 cpl_vector *res = NULL;
1312 cpl_ensure(v != NULL, CPL_ERROR_NULL_INPUT, res);
1313 int size = cpl_vector_get_size(v);
1314
1315 int nbnan = 0;
1316 for (int i = 0; i < size; i++) {
1317 double val = cpl_vector_get(v, i);
1318 if (isnan(val)) {
1319 nbnan++;
1320 }
1321 }
1322 int nsize = size - nbnan;
1323 if (nsize > 0) {
1324 res = cpl_vector_new(nsize);
1325 int j = 0;
1326 for (int i = 0; i < size; i++) {
1327 double val = cpl_vector_get(v, i);
1328 if (!isnan(val)) {
1329 cpl_vector_set(res, j, val);
1330 j++;
1331 }
1332 }
1333 }
1334 return res;
1335}
1336/*----------------------------------------------------------------------------*/
1342cpl_bivector *
1344{
1345 cpl_bivector *res = NULL;
1346 cpl_ensure(v != NULL, CPL_ERROR_NULL_INPUT, NULL);
1347
1348 int size = cpl_bivector_get_size(v);
1349 cpl_vector *y = cpl_bivector_get_y(v);
1350 cpl_vector *x = cpl_bivector_get_x(v);
1351 int nbnan = 0;
1352 for (int i = 0; i < size; i++) {
1353 double xval = cpl_vector_get(x, i);
1354 double yval = cpl_vector_get(y, i);
1355 if (isnan(xval) || isnan(yval)) {
1356 nbnan++;
1357 }
1358 }
1359 int nsize = size - nbnan;
1360 if (nsize > 0) {
1361 res = cpl_bivector_new(nsize);
1362 cpl_vector *ry = cpl_bivector_get_y(res);
1363 cpl_vector *rx = cpl_bivector_get_x(res);
1364
1365 int j = 0;
1366 for (int i = 0; i < size; i++) {
1367 double xval = cpl_vector_get(x, i);
1368 double yval = cpl_vector_get(y, i);
1369 if (!isnan(xval) && !isnan(yval)) {
1370 cpl_vector_set(ry, j, xval);
1371 cpl_vector_set(rx, j, yval);
1372 j++;
1373 }
1374 }
1375 }
1376 return res;
1377}
1378
1379/*----------------------------------------------------------------------------*/
1385hdrl_image *
1387{
1388 cpl_ensure(image != NULL, CPL_ERROR_NULL_INPUT, NULL);
1389
1390 hdrl_image *res = NULL;
1391
1392 int nx = hdrl_image_get_size_x(image);
1393 int ny = hdrl_image_get_size_y(image);
1394
1395 res = hdrl_image_new(nx, 1);
1396
1397 for (int i = 1; i <= nx; i++) {
1398 hdrl_image *col = hdrl_image_extract(image, i, 1, i, ny);
1399
1400 hdrl_value median = hdrl_image_get_median(col);
1401 hdrl_image_set_pixel(res, i, 1, median);
1402 hdrl_image_delete(col);
1403 }
1404 return res;
1405}
1406
1407/*----------------------------------------------------------------------------*/
1416double
1417moo_sky_distance(double alpha1, double delta1, double alpha2, double delta2)
1418{
1419 double res = 0;
1420
1421 double alpha1_rad = CPL_MATH_RAD_DEG * alpha1;
1422 double delta1_rad = CPL_MATH_RAD_DEG * delta1;
1423 double alpha2_rad = CPL_MATH_RAD_DEG * alpha2;
1424 double delta2_rad = CPL_MATH_RAD_DEG * delta2;
1425
1426 res =
1427 acos(sin(delta1_rad) * sin(delta2_rad) +
1428 cos(delta1_rad) * cos(delta2_rad) * cos(alpha1_rad - alpha2_rad));
1429 return res;
1430}
1431
1432moo_spline *
1433moo_spline_create(cpl_bivector *data)
1434{
1435 cpl_ensure(data != NULL, CPL_ERROR_NULL_INPUT, NULL);
1436
1437 int size = cpl_bivector_get_size(data);
1438 double *x = cpl_bivector_get_x_data(data);
1439 double *y = cpl_bivector_get_y_data(data);
1440
1441 moo_spline *res = cpl_calloc(1, sizeof(moo_spline));
1442 res->min = x[0];
1443 res->max = x[size - 1];
1444 res->acc = gsl_interp_accel_alloc();
1445 res->spline = gsl_spline_alloc(gsl_interp_cspline, size);
1446 gsl_spline_init(res->spline, x, y, size);
1447
1448 return res;
1449}
1450double
1451moo_spline_eval(moo_spline *self, double x)
1452{
1453 double res = NAN;
1454 cpl_ensure(self != NULL, CPL_ERROR_NULL_INPUT, 0.0);
1455 if (x >= self->min && x <= self->max) {
1456 res = gsl_spline_eval(self->spline, x, self->acc);
1457 }
1458 return res;
1459}
1460
1461void
1462moo_spline_delete(moo_spline *self)
1463{
1464 if (self != NULL) {
1465 gsl_interp_accel_free(self->acc);
1466 gsl_spline_free(self->spline);
1467 cpl_free(self);
1468 }
1469}
1470
1471
1472static int
1473_moo_gsl_fit_wmul(const double *x,
1474 const double *w,
1475 const double *y,
1476 const double *ey,
1477 const size_t n,
1478 double *c1,
1479 double *sig_c1,
1480 double *cov_11,
1481 double *chisq)
1482{
1483 /* compute the weighted means and weighted deviations from the means */
1484
1485 /* wm denotes a "weighted mean", wm(f) = (sum_i w_i f_i) / (sum_i w_i) */
1486
1487 double W = 0, wm_x = 0, wm_y = 0, wm_dx2 = 0, wm_dxdy = 0;
1488 double sig_wm_y2 = 0, sig_wm_dxdy2 = 0; /* FRED DEBUG */
1489
1490 /* ey[i] denotes l'erreur sur y - FRED DEBUG */
1491
1492 size_t i;
1493
1494 for (i = 0; i < n; i++) {
1495 const double wi = w[i];
1496
1497 if (wi > 0) {
1498 W += wi;
1499 wm_x += (x[i] - wm_x) * (wi / W);
1500 wm_y += (y[i] - wm_y) * (wi / W);
1501 sig_wm_y2 += (wi * wi * ey[i] * ey[i]); /* FRED DEBUG */
1502 }
1503 }
1504
1505 W = 0; /* reset the total weight */
1506
1507 for (i = 0; i < n; i++) {
1508 const double wi = w[i];
1509
1510 if (wi > 0) {
1511 const double dx = x[i] - wm_x;
1512 const double dy = y[i] - wm_y;
1513
1514 W += wi;
1515 wm_dx2 += (dx * dx - wm_dx2) * (wi / W);
1516 wm_dxdy += (dx * dy - wm_dxdy) * (wi / W);
1517 sig_wm_dxdy2 += wi * wi * dx * dx *
1518 (ey[i] * ey[i] + sig_wm_y2); /* FRED DEBUG */
1519 }
1520 }
1521
1522 /* In terms of y = b x */
1523
1524 {
1525 double d2 = 0;
1526 double b = (wm_x * wm_y + wm_dxdy) / (wm_x * wm_x + wm_dx2);
1527
1528 double sig_b = sqrt(wm_x * wm_x * sig_wm_y2 + sig_wm_dxdy2) /
1529 (wm_x * wm_x + wm_dx2); /* FRED DEBUG */
1530
1531 *c1 = b;
1532 *sig_c1 = sig_b;
1533 *cov_11 = 1 / (W * (wm_x * wm_x + wm_dx2));
1534
1535 /* Compute chi^2 = \sum w_i (y_i - b * x_i)^2 */
1536
1537 for (i = 0; i < n; i++) {
1538 const double wi = w[i];
1539
1540 if (wi > 0) {
1541 const double dx = x[i] - wm_x;
1542 const double dy = y[i] - wm_y;
1543 const double d = (wm_y - b * wm_x) + (dy - b * dx);
1544 d2 += wi * d * d;
1545 }
1546 }
1547
1548 *chisq = d2;
1549 }
1550
1551 return GSL_SUCCESS;
1552}
1553
1554/*----------------------------------------------------------------------------*/
1566cpl_error_code
1567moo_fit_mul(const cpl_vector *vx,
1568 const cpl_vector *vw,
1569 const cpl_vector *vy,
1570 const cpl_vector *vy_err,
1571 double *c,
1572 double *sig_c)
1573{
1574 cpl_ensure_code(vx != NULL, CPL_ERROR_NULL_INPUT);
1575 cpl_ensure_code(vy != NULL, CPL_ERROR_NULL_INPUT);
1576 cpl_ensure_code(vw != NULL, CPL_ERROR_NULL_INPUT);
1577
1578 int n = cpl_vector_get_size(vx);
1579
1580 const double *x = cpl_vector_get_data_const(vx);
1581 const double *y = cpl_vector_get_data_const(vy);
1582 const double *yerr = cpl_vector_get_data_const(vy_err);
1583 const double *w = cpl_vector_get_data_const(vw);
1584
1585 double cov11, chisq;
1586
1587 _moo_gsl_fit_wmul(x, w, y, yerr, n, c, sig_c, &cov11, &chisq);
1588
1589
1590 return CPL_ERROR_NONE;
1591}
1592
1593static cpl_polynomial *
1594_moo_fit(cpl_vector *v, int polyorder)
1595{
1596 gsl_matrix *X = NULL, *cov = NULL;
1597 gsl_vector *y = NULL, *c = NULL;
1598 gsl_multifit_linear_workspace *work = NULL;
1599 cpl_size pows[1];
1600
1601 cpl_ensure(v != NULL, CPL_ERROR_NULL_INPUT, NULL);
1602
1603 cpl_polynomial *pol = cpl_polynomial_new(1);
1604 int i, j;
1605 double chisq;
1606 double *data = cpl_vector_get_data(v);
1607 int size = cpl_vector_get_size(v);
1608
1609 int nbcoeffs = polyorder + 1;
1610 X = gsl_matrix_alloc(size, nbcoeffs);
1611 y = gsl_vector_alloc(size);
1612 c = gsl_vector_alloc(nbcoeffs);
1613 cov = gsl_matrix_alloc(nbcoeffs, nbcoeffs);
1614
1615 int nbnan = 0;
1616 for (i = 0; i < size; i++) {
1617 double d = data[i];
1618 if (!isnan(d)) {
1619 for (j = 0; j < polyorder + 1; j++) {
1620 gsl_matrix_set(X, i, j, pow(i, j));
1621 }
1622 gsl_vector_set(y, i, d);
1623 }
1624 else {
1625 nbnan++;
1626 }
1627 }
1628
1629 if ((size - nbnan) > (nbcoeffs + 1)) {
1630 work = gsl_multifit_linear_alloc(size, nbcoeffs);
1631 gsl_multifit_linear(X, y, c, cov, &chisq, work);
1632
1633 for (j = 0; j < polyorder + 1; j++) {
1634 pows[0] = j;
1635 cpl_polynomial_set_coeff(pol, pows, gsl_vector_get(c, j));
1636 }
1637 gsl_multifit_linear_free(work);
1638 }
1639 gsl_matrix_free(X);
1640 gsl_vector_free(y);
1641 gsl_vector_free(c);
1642 gsl_matrix_free(cov);
1643 return pol;
1644}
1645/*----------------------------------------------------------------------------*/
1655cpl_vector *
1656moo_savgol_filter(cpl_vector *v, int window_length, int polyorder)
1657{
1658 cpl_ensure(v != NULL, CPL_ERROR_NULL_INPUT, NULL);
1659 cpl_ensure(polyorder >= 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
1660 cpl_ensure(polyorder < window_length, CPL_ERROR_ILLEGAL_INPUT, NULL);
1661 cpl_ensure(window_length >= 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
1662
1663 int size = cpl_vector_get_size(v);
1664 cpl_vector *result = cpl_vector_new(size);
1665 double *resdata = cpl_vector_get_data(result);
1666
1667 cpl_polynomial *fit = NULL;
1668 cpl_vector *data = NULL;
1669 cpl_size pows[1];
1670
1671 int hsize = (window_length - 1) / 2;
1672 /* first values */
1673 data = cpl_vector_extract(v, 0, window_length - 1, 1);
1674
1675 fit = _moo_fit(data, polyorder);
1676
1677 for (int i = 0; i < hsize; i++) {
1678 double yval = 0;
1679 for (int j = 0; j < polyorder + 1; j++) {
1680 pows[0] = j;
1681 double coef = cpl_polynomial_get_coeff(fit, pows);
1682 yval += coef * pow(i, j);
1683 }
1684 resdata[i] = yval;
1685 }
1686 cpl_vector_delete(data);
1687 cpl_polynomial_delete(fit);
1688 /* end values */
1689 data = cpl_vector_extract(v, size - window_length, size - 1, 1);
1690 fit = _moo_fit(data, polyorder);
1691
1692 for (int i = hsize; i < window_length; i++) {
1693 double yval = 0;
1694 for (int j = 0; j < polyorder + 1; j++) {
1695 pows[0] = j;
1696 double coef = cpl_polynomial_get_coeff(fit, pows);
1697 yval += coef * pow(i, j);
1698 }
1699 resdata[size - window_length + i] = yval;
1700 }
1701 cpl_vector_delete(data);
1702 cpl_polynomial_delete(fit);
1703
1704 /* loop */
1705 for (int i = hsize; i < size - hsize; i++) {
1706 data =
1707 cpl_vector_extract(v, i - hsize, i - hsize - 1 + window_length, 1);
1708 fit = _moo_fit(data, polyorder);
1709 double yval = 0;
1710 for (int j = 0; j < polyorder + 1; j++) {
1711 pows[0] = j;
1712 double coef = cpl_polynomial_get_coeff(fit, pows);
1713 yval += coef * pow(hsize, j);
1714 }
1715 cpl_vector_delete(data);
1716 cpl_polynomial_delete(fit);
1717 resdata[i] = yval;
1718 }
1719 return result;
1720}
1721
1722/*----------------------------------------------------------------------------*/
1729cpl_vector *
1730moo_median_filter(cpl_vector *v, int winhsize)
1731{
1732 cpl_ensure(v != NULL, CPL_ERROR_NULL_INPUT, NULL);
1733 cpl_ensure(winhsize >= 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
1734
1735 cpl_vector *result = NULL;
1736 int size = cpl_vector_get_size(v);
1737 int K = winhsize * 2 + 1;
1738
1739 result = cpl_vector_new(size);
1740
1741 double *windata = cpl_calloc(K, sizeof(double));
1742 /* generate input signal */
1743 for (int i = 0; i < size; ++i) {
1744 int idx = 0;
1745
1746 int kmin = i - winhsize;
1747 if (kmin < 0) {
1748 kmin = 0;
1749 }
1750
1751 int kmax = i + winhsize;
1752 if (kmax >= size) {
1753 kmax = size - 1;
1754 }
1755 for (int k = kmin; k <= kmax; k++) {
1756 double val = cpl_vector_get(v, k);
1757 if (!isnan(val)) {
1758 windata[idx] = val;
1759 idx++;
1760 }
1761 }
1762 if (idx > 0) {
1763 cpl_vector *winv = cpl_vector_wrap(idx, windata);
1764
1765 double med = cpl_vector_get_median(winv);
1766 cpl_vector_set(result, i, med);
1767 cpl_vector_unwrap(winv);
1768 }
1769 else {
1770 cpl_vector_set(result, i, NAN);
1771 }
1772 }
1773
1774 cpl_free(windata);
1775
1776 return result;
1777}
1778
1779/*----------------------------------------------------------------------------*/
1790cpl_image *
1791moo_imagelist_collapse_bitwiseor(cpl_imagelist *list, hdrl_imagelist *datalist)
1792{
1793 cpl_image *result = NULL;
1794 int size = 0;
1795
1796 cpl_errorstate prestate = cpl_errorstate_get();
1797
1798 cpl_ensure(list != NULL, CPL_ERROR_NULL_INPUT, NULL);
1799 cpl_ensure(datalist != NULL, CPL_ERROR_NULL_INPUT, NULL);
1800 size = cpl_imagelist_get_size(list);
1801 cpl_ensure(size > 0, CPL_ERROR_NULL_INPUT, NULL);
1802
1803 cpl_image *first = cpl_imagelist_get(list, 0);
1804 int nx = cpl_image_get_size_x(first);
1805 int ny = cpl_image_get_size_y(first);
1806
1807 result = cpl_image_new(nx, ny, CPL_TYPE_INT);
1808
1809 for (int y = 1; y <= ny; y++) {
1810 for (int x = 1; x <= nx; x++) {
1811 int v = 0;
1812 int av = 0;
1813 int nbcomb = 0;
1814 for (int k = 0; k < size; k++) {
1815 cpl_image *qual = cpl_imagelist_get(list, k);
1816 hdrl_image *data = hdrl_imagelist_get(datalist, k);
1817 int rej;
1818 int qv = cpl_image_get(qual, x, y, &rej);
1819 hdrl_image_get_pixel(data, x, y, &rej);
1820 av |= qv;
1821 if (rej == 0) {
1822 v |= qv;
1823 nbcomb++;
1824 }
1825 }
1826 if (nbcomb == 0) {
1827 v = av;
1828 }
1829 cpl_image_set(result, x, y, v);
1830 }
1831 }
1832
1833
1834 if (!cpl_errorstate_is_equal(prestate)) {
1835 cpl_image_delete(result);
1836 result = NULL;
1837 cpl_errorstate_dump(prestate, CPL_FALSE, cpl_errorstate_dump_one);
1838 // Recover from the error(s) (Reset to prestate))
1839 cpl_errorstate_set(prestate);
1840 }
1841 return result;
1842}
1843
1844/*----------------------------------------------------------------------------*/
1862double
1863moo_image_get_ron(cpl_image *diff,
1864 int llx,
1865 int lly,
1866 int urx,
1867 int ury,
1868 int nb_boxes,
1869 int box_hsize,
1870 double max_error_frac,
1871 int max_niter)
1872{
1873 double ron = 0.0;
1874 double ronerr = 0.0;
1875 int niter = 0;
1876
1877 cpl_ensure(diff, CPL_ERROR_NULL_INPUT, 0.0);
1878 cpl_size *window = (cpl_size *)cpl_calloc(sizeof(cpl_size), 4);
1879 window[0] += llx;
1880 window[1] += urx;
1881 window[2] += lly;
1882 window[3] += ury;
1883
1884 do {
1885 cpl_flux_get_noise_window(diff, window, box_hsize, nb_boxes, &ron,
1886 &ronerr);
1887 niter++;
1888 } while (ronerr > max_error_frac * ron && niter < max_niter);
1889
1890 cpl_free(window);
1891 return ron;
1892}
1893
1894/*----------------------------------------------------------------------------*/
1905cpl_error_code
1906moo_image_get_quartile(cpl_image *image, double *qmin, double *qmax)
1907{
1908 cpl_ensure_code(image != NULL, CPL_ERROR_NULL_INPUT);
1909 cpl_ensure_code(qmin != NULL, CPL_ERROR_NULL_INPUT);
1910 cpl_ensure_code(qmax != NULL, CPL_ERROR_NULL_INPUT);
1911
1912 int nx = cpl_image_get_size_x(image);
1913 int ny = cpl_image_get_size_y(image);
1914 int size = nx * ny;
1915 cpl_image *stats = cpl_image_cast(image, CPL_TYPE_DOUBLE);
1916 cpl_image_reject_value(stats, CPL_VALUE_NAN);
1917 int nb = cpl_image_count_rejected(stats);
1918 int nsize = size - nb;
1919 cpl_vector *statsv = cpl_vector_new(nsize);
1920 double *data = cpl_image_get_data_double(stats);
1921 int nbnan = 0;
1922 for (int i = 0; i < size; i++) {
1923 double v = data[i];
1924 if (isnan(v)) {
1925 nbnan++;
1926 }
1927 else {
1928 cpl_vector_set(statsv, i - nbnan, v);
1929 }
1930 }
1931 double *stats_data = cpl_vector_get_data(statsv);
1932 gsl_sort(stats_data, 1, nsize);
1933 double upperq =
1934 gsl_stats_quantile_from_sorted_data(stats_data, 1, nsize, 0.8415);
1935 double lowerq =
1936 gsl_stats_quantile_from_sorted_data(stats_data, 1, nsize, 0.1585);
1937 *qmin = lowerq;
1938 *qmax = upperq;
1939
1940 cpl_vector_delete(statsv);
1941 cpl_image_delete(stats);
1942 return CPL_ERROR_NONE;
1943}
1944
1945/* ---------------------------------------------------------------------------*/
1946#define MACHEP DBL_EPSILON
1947#define MAXLOG log(FLT_MAX)
1948#define mtherr(x, y)
1949
1950static double igam(double, double);
1951static double igamc(double, double);
1952
1953static double big = 4.503599627370496e15;
1954static double biginv = 2.22044604925031308085e-16;
1955
1956static double
1957igam(double a, double x)
1958{
1959 double ans, ax, c, r;
1960
1961 /* Check zero integration limit first */
1962 if (x == 0) {
1963 return (0.);
1964 }
1965
1966 if ((x < 0) || (a <= 0)) {
1967 mtherr("igam", DOMAIN);
1968 return (NAN);
1969 }
1970
1971 if ((x > 1.) && (x > a)) {
1972 return (1. - igamc(a, x));
1973 }
1974
1975 /* Compute x**a * exp(-x) / gamma(a) */
1976 ax = a * log(x) - x - lgamma(a);
1977 if (ax < -MAXLOG) {
1978 mtherr("igam", UNDERFLOW);
1979 return (0.);
1980 }
1981 ax = exp(ax);
1982
1983 /* power series */
1984 r = a;
1985 c = 1.;
1986 ans = 1.;
1987
1988 do {
1989 r += 1.;
1990 c *= x / r;
1991 ans += c;
1992
1993 } while (c / ans > MACHEP);
1994
1995 return (ans * ax / a);
1996}
1997
1998static double
1999igamc(double a, double x)
2000{
2001 double ans, ax, c, r, t, y, z;
2002 double pkm1, pkm2, qkm1, qkm2;
2003
2004 if ((x < 0) || (a <= 0)) {
2005 mtherr("igamc", DOMAIN);
2006 return (NAN);
2007 }
2008
2009 if ((x < 1.0) || (x < a)) {
2010 return (1. - igam(a, x));
2011 }
2012 ax = a * log(x) - x - lgamma(a);
2013
2014 if (ax < -MAXLOG) {
2015 mtherr("igamc", UNDERFLOW);
2016 return (0.);
2017 }
2018 ax = exp(ax);
2019
2020 /* continued fraction */
2021 y = 1. - a;
2022 z = x + y + 1.;
2023 c = 0.;
2024 pkm2 = 1.;
2025 qkm2 = x;
2026 pkm1 = x + 1.;
2027 qkm1 = z * x;
2028 ans = pkm1 / qkm1;
2029
2030 do {
2031 c += 1.;
2032 y += 1.;
2033 z += 2.;
2034
2035 double yc = y * c;
2036 double pk = pkm1 * z - pkm2 * yc;
2037 double qk = qkm1 * z - qkm2 * yc;
2038
2039 if (qk != 0) {
2040 r = pk / qk;
2041 t = fabs((ans - r) / r);
2042 ans = r;
2043 }
2044 else {
2045 t = 1.;
2046 }
2047
2048 pkm2 = pkm1;
2049 pkm1 = pk;
2050 qkm2 = qkm1;
2051 qkm1 = qk;
2052
2053 if (fabs(pk) > big) {
2054 pkm2 *= biginv;
2055 pkm1 *= biginv;
2056 qkm2 *= biginv;
2057 qkm1 *= biginv;
2058 }
2059
2060 } while (t > MACHEP);
2061
2062 return (ans * ax);
2063}
2064
2065static cpl_image *
2066bpm_from_rel(cpl_image *img, double kappa_low, double kappa_high, int mad)
2067{
2068 cpl_image *r;
2069 double m, s;
2070 if (mad) {
2071 m = cpl_image_get_mad(img, &s);
2072 s *= CPL_MATH_STD_MAD;
2073 s = CX_MAX(DBL_EPSILON, s);
2074 }
2075 else {
2076 m = cpl_image_get_mean(img);
2077 s = cpl_image_get_stdev(img);
2078 }
2079 cpl_mask *bpm = cpl_mask_threshold_image_create(img, m - s * kappa_low,
2080 m + s * kappa_high);
2081 cpl_mask_not(bpm);
2082 r = cpl_image_new_from_mask(bpm);
2083 cpl_mask_delete(bpm);
2084 return r;
2085}
2086cpl_error_code
2087moo_hdrl_bpm_fit_compute(const hdrl_parameter *par,
2088 const hdrl_imagelist *data,
2089 const cpl_vector *sample_pos,
2090 cpl_image **out_mask)
2091{
2092 cpl_error_code status = CPL_ERROR_NONE;
2093 cpl_image *out_chi2 = NULL, *out_dof = NULL;
2094 hdrl_imagelist *out_coef = NULL;
2095
2096 status = hdrl_bpm_fit_parameter_verify(par);
2097 if (status != CPL_ERROR_NONE)
2098 return status;
2099
2100 int degree = hdrl_bpm_fit_parameter_get_degree(par);
2101
2102 /* TODO check for 0 error,
2103 * in that case set to 1 and do not allow pval */
2104
2105 status = hdrl_fit_polynomial_imagelist(data, sample_pos, degree, &out_coef,
2106 &out_chi2, &out_dof);
2107
2108 if (status) {
2109 return cpl_error_set_message(cpl_func, CPL_ERROR_FILE_NOT_FOUND,
2110 "Fit failed");
2111 }
2112
2113 if (cpl_image_count_rejected(out_chi2) ==
2114 (cpl_image_get_size_x(out_chi2) * cpl_image_get_size_y(out_chi2))) {
2115 cpl_msg_error(cpl_func,
2116 "Too few good pixels to fit polynomial of "
2117 "degree %d in all pixels",
2118 degree);
2119
2120 goto end;
2121 }
2122
2123 {
2124 cpl_image *bpm = NULL;
2125 double pval = hdrl_bpm_fit_parameter_get_pval(par);
2126 double rel_chi_l = hdrl_bpm_fit_parameter_get_rel_chi_low(par);
2127 double rel_chi_h = hdrl_bpm_fit_parameter_get_rel_chi_high(par);
2128 double rel_coef_l = hdrl_bpm_fit_parameter_get_rel_coef_low(par);
2129 double rel_coef_h = hdrl_bpm_fit_parameter_get_rel_coef_high(par);
2130
2131 if (rel_chi_l >= 0) {
2132 /* chi is symmetric */
2133 cpl_image_power(out_chi2, 0.5);
2134 bpm = bpm_from_rel(out_chi2, rel_chi_l, rel_chi_h, 1);
2135 }
2136 else if (rel_coef_l >= 0) {
2137 for (cpl_size i = 0; i < hdrl_imagelist_get_size(out_coef); i++) {
2138 hdrl_image *coef = hdrl_imagelist_get(out_coef, i);
2139 cpl_image *b = bpm_from_rel(hdrl_image_get_image(coef),
2140 rel_coef_l, rel_coef_h, 1);
2141
2142 /* bits of bpm defines which coefficient "bad" */
2143 if (bpm) {
2144 cpl_image_multiply_scalar(b, pow(2, i));
2145 cpl_image_add(bpm, b);
2146 cpl_image_delete(b);
2147 }
2148 else {
2149 bpm = b;
2150 }
2151 }
2152 }
2153 else if (pval >= 0) {
2154 pval /= 100.;
2155 bpm = cpl_image_new(cpl_image_get_size_x(out_chi2),
2156 cpl_image_get_size_y(out_chi2), CPL_TYPE_INT);
2157 int *md = cpl_image_get_data_int(bpm);
2158 hdrl_data_t *cd = cpl_image_get_data(out_chi2);
2159 hdrl_data_t *dd = cpl_image_get_data(out_dof);
2160 for (size_t i = 0; i < (size_t)cpl_image_get_size_x(out_chi2) *
2161 cpl_image_get_size_y(out_chi2);
2162 i++) {
2163 double pv = igamc(dd[i] / 2., cd[i] / 2.);
2164 md[i] = (pv < pval);
2165 }
2166 }
2167 *out_mask = bpm;
2168 }
2169
2170end:
2171 hdrl_imagelist_delete(out_coef);
2172 cpl_image_delete(out_chi2);
2173 cpl_image_delete(out_dof);
2174
2175 return cpl_error_get_code();
2176}
2177/*----------------------------------------------------------------------------*/
2188double
2189moo_vector_get_dersnr(const cpl_vector *ve)
2190{
2191 double snr = 0.0;
2192 int size = cpl_vector_get_size(ve);
2193
2194 if (size > 4) {
2195 cpl_vector *sve = cpl_vector_duplicate(ve);
2196 double signal = cpl_vector_get_median(sve);
2197 cpl_vector_delete(sve);
2198 cpl_vector *ve1 = cpl_vector_extract(ve, 2, size - 3, 1);
2199 cpl_vector *ve2 = cpl_vector_extract(ve, 0, size - 5, 1);
2200 cpl_vector *ve3 = cpl_vector_extract(ve, 4, size - 1, 1);
2201 int size1 = cpl_vector_get_size(ve1);
2202
2203 for (int i = 0; i < size1; i++) {
2204 double v1 = cpl_vector_get(ve1, i);
2205 double v2 = cpl_vector_get(ve2, i);
2206 double v3 = cpl_vector_get(ve3, i);
2207 double v = fabs(2 * v1 - v2 - v3);
2208 cpl_vector_set(ve1, i, v);
2209 }
2210 double noise = 0.6052697 * cpl_vector_get_median(ve1);
2211
2212 snr = signal / noise;
2213
2214 cpl_vector_delete(ve1);
2215 cpl_vector_delete(ve2);
2216 cpl_vector_delete(ve3);
2217 }
2218
2219 return snr;
2220}
2221/*----------------------------------------------------------------------------*/
2232int
2233moo_string_is_strictly_equal(const char *a, const char *b)
2234{
2235 int res = 0;
2236 cpl_ensure(a != NULL, CPL_ERROR_NULL_INPUT, 0);
2237 cpl_ensure(b != NULL, CPL_ERROR_NULL_INPUT, 0);
2238
2239 int sizea = strlen(a);
2240 int sizeb = strlen(b);
2241
2242 if (sizea != sizeb) {
2243 res = 0;
2244 }
2245 else {
2246 res = 1;
2247 for (int i = 0; i < sizea; i++) {
2248 if (a[i] != b[i]) {
2249 res = 0;
2250 break;
2251 }
2252 }
2253 }
2254 return res;
2255}
2256
2257moo_date
2258moo_get_date_from_string(const char *string)
2259{
2260 moo_date date;
2261 int year, month, day, hour, min;
2262 float sec;
2263 int nb = sscanf(string, "%d-%d-%dT%d:%d:%f", &year, &month, &day, &hour,
2264 &min, &sec);
2265 date.day = day;
2266 date.year = year;
2267 date.month = month;
2268 date.hour = hour;
2269 date.min = min;
2270 date.sec = sec;
2271 return date;
2272}
2273/*----------------------------------------------------------------------------*/
2284cpl_size
2286 const char *name,
2287 const char *string)
2288{
2289 int nrow = cpl_table_get_nrow(table);
2290 cpl_size nbsel = 0;
2291 const char **names = cpl_table_get_data_string_const(table, name);
2292 for (int i = 0; i < nrow; i++) {
2293 if (!cpl_table_is_selected(table, i)) {
2294 if (moo_string_is_strictly_equal(names[i], string)) {
2295 cpl_table_select_row(table, i);
2296 nbsel++;
2297 }
2298 }
2299 else {
2300 nbsel++;
2301 }
2302 }
2303 return nbsel;
2304}
2305/*----------------------------------------------------------------------------*/
2316cpl_size
2318 const char *name,
2319 const char *string)
2320{
2321 char *regexp = cpl_sprintf("^%s$", string);
2322 int size = cpl_table_and_selected_string(table, name, CPL_EQUAL_TO, regexp);
2323 cpl_free(regexp);
2324 return size;
2325}
2326/*******************************************************************************/
2327/* only for gsl2.6 */
2328#if 0
2329cpl_vector* moo_median_filter(cpl_vector* v,int winhsize)
2330{
2331 cpl_ensure(v!=NULL,CPL_ERROR_NULL_INPUT,NULL);
2332 cpl_ensure(winhsize>=0,CPL_ERROR_ILLEGAL_INPUT,NULL);
2333
2334 cpl_vector* result = NULL;
2335 int size = cpl_vector_get_size(v);
2336 int K = winhsize*2+1;
2337
2338
2339 gsl_filter_median_workspace *median_p = gsl_filter_median_alloc(K);
2340 gsl_vector *x = gsl_vector_alloc(size); /* input vector */
2341 gsl_vector *y_median = gsl_vector_alloc(size); /* median filtered output */
2342
2343 /* generate input signal */
2344 for (int i = 0; i < size; ++i)
2345 {
2346 double val = cpl_vector_get(v,i);
2347 gsl_vector_set(x, i, val);
2348 }
2349
2350 gsl_filter_median(GSL_FILTER_END_PADZERO, x, y_median, median_p);
2351
2352 result = cpl_vector_new(size);
2353
2354 for (int i = 0; i < size; ++i)
2355 {
2356 double val = gsl_vector_get(y_median,i);
2357 cpl_vector_set(result, i, val);
2358 }
2359
2360 gsl_vector_free(x);
2361 gsl_vector_free(y_median);
2362 gsl_filter_median_free(median_p);
2363
2364 return result;
2365}
2366
2367
2368#endif
2369
#define MOO_BADPIX_GOOD
Definition: moo_badpix.h:36
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_sky_distance(double alpha1, double delta1, double alpha2, double delta2)
Compute sky distance (in rad)
Definition: moo_utils.c:1417
cpl_error_code moo_gaussian_eval_inv(double y, double x0, double sigma, double offset, double area, double *x1, double *x2)
Find the x positions of the gaussian at the given y position.
Definition: moo_utils.c:373
cpl_size moo_table_or_selected_sequal_string(cpl_table *table, const char *name, const char *string)
Select from unselected table rows, by comparing column values with a constant.
Definition: moo_utils.c:2285
hdrl_image * moo_image_collapse_median_create(hdrl_image *image)
Collapse row of an image using a median and compute associate error.
Definition: moo_utils.c:1386
double moo_vector_get_min(const cpl_vector *v, int *flags)
Find minimum values in a vector using flags.
Definition: moo_utils.c:979
cpl_bivector * moo_bivector_filter_nan(cpl_bivector *v)
Create new bi vector with nan values filter.
Definition: moo_utils.c:1343
double moo_gaussian_eval(double x, double x0, double sigma, double offset, double area)
Evaluate the gaussian at the given x position.
Definition: moo_utils.c:348
moo_tcheby_polynomial * moo_tcheby_polynomial_fit(cpl_bivector *data, int degree, double xmin, double xmax)
Computes Tchebitchev transformation of data.
Definition: moo_utils.c:573
cpl_mask * moo_kappa_sigma_clipping(cpl_image *sigma_img, int niter, double kappa, double cdiff, double maxfrac)
Compute mask of rejected pixels using kappa sigma algorithm.
Definition: moo_utils.c:163
cpl_error_code moo_gaussian_fit(cpl_bivector *points, cpl_fit_mode fit_pars, double *center, double *width, double *background, double *area)
Fit the data with a gaussian.
Definition: moo_utils.c:309
cpl_vector * moo_vector_filter_nan(cpl_vector *v)
Create new vector with nan values filter.
Definition: moo_utils.c:1309
cpl_image * moo_imagelist_collapse_bitwiseor(cpl_imagelist *list, hdrl_imagelist *datalist)
get the QUAL resulting in a bitwise OR operation on the QUAL list
Definition: moo_utils.c:1791
cpl_size moo_table_and_selected_sequal_string(cpl_table *table, const char *name, const char *string)
Select from unselected table rows, by comparing column values with a constant.
Definition: moo_utils.c:2317
cpl_error_code moo_find_threshold_limits(cpl_bivector *points, double threshold, double *xmin, double *xmax)
Find threshold limits of a 1D signal.
Definition: moo_utils.c:423
moo_tcheby2d_polynomial * moo_tcheby2d_polynomial_fit(cpl_vector *in_x, int xdegree, double xmin, double xmax, cpl_vector *in_y, int ydegree, double ymin, double ymax, cpl_vector *in_l, double lmin, double lmax)
Computes Tchebitchev transformation of data.
Definition: moo_utils.c:667
double moo_image_get_ron(cpl_image *diff, int llx, int lly, int urx, int ury, int nb_boxes, int box_hsize, double max_error_frac, int max_niter)
compute ron in a diff image using boxes
Definition: moo_utils.c:1863
double moo_vector_get_max(const cpl_vector *v, int *flags)
Find maximum values in a vector using flags.
Definition: moo_utils.c:1003
const char * moo_get_license(void)
Get the pipeline copyright and license.
Definition: moo_utils.c:81
cpl_error_code moo_fit_mul(const cpl_vector *vx, const cpl_vector *vw, const cpl_vector *vy, const cpl_vector *vy_err, double *c, double *sig_c)
This function computes the best-fit linear regression coefficient c1 of the model Y = c_1 X for the w...
Definition: moo_utils.c:1567
cpl_error_code moo_interpolate_linear(cpl_bivector *fout, cpl_vector *fout_errs, int *fout_qual, const cpl_bivector *fref, const cpl_vector *fref_errs, const int *fref_qual)
Linear interpolation of a 1d-function.
Definition: moo_utils.c:1160
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_dersnr(const cpl_vector *ve)
This function computes the signal to noise ratio DER_SNR following the definition set forth by the Sp...
Definition: moo_utils.c:2189
cpl_error_code moo_image_get_quartile(cpl_image *image, double *qmin, double *qmax)
compute first and last quartile from an image
Definition: moo_utils.c:1906
double moo_vector_get_percentile(cpl_vector *v, double f)
Get percentile of input vector.
Definition: moo_utils.c:1263
cpl_error_code moo_barycenter_fit(cpl_bivector *points, double *center, double *width)
Fit positions using weighted fluxes.
Definition: moo_utils.c:254
int moo_string_is_strictly_equal(const char *a, const char *b)
This function compares to string to see if the two string are stricly equal.
Definition: moo_utils.c:2233
cpl_image * moo_compute_sigma_map(cpl_imagelist *list, cpl_imagelist *qlist, cpl_image *img)
Compute image of sigma variation.
Definition: moo_utils.c:100