CR2RE Pipeline Reference Manual 1.6.10
hdrl_fringe.c
1/*
2 * This file is part of the HDRL
3 * Copyright (C) 2015 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 St, Fifth Floor, Boston, MA 02110-1301 USA
18 *
19 * hdrl_fringe.c
20 *
21 * Created on: May 11, 2015
22 * Author: agabasch
23 */
24
25#ifdef HAVE_CONFIG_H
26#include <config.h>
27#endif
28
29/*-----------------------------------------------------------------------------
30 Includes
31 -----------------------------------------------------------------------------*/
32
33#include <hdrl_fringe.h>
34#include <hdrl_prototyping.h>
35#include <math.h>
36
37
38
39
40/*----------------------------------------------------------------------------*/
114/*----------------------------------------------------------------------------*/
115
118/*---------------------------------------------------------------------------*/
158/*---------------------------------------------------------------------------*/
159cpl_error_code
160hdrl_fringe_compute(hdrl_imagelist* ilist_fringe, const cpl_imagelist * ilist_obj,
161 const cpl_mask* stat_mask, const hdrl_parameter* collapse_params,
162 hdrl_image** master, cpl_image** contrib_map,
163 cpl_table ** qctable)
164{
165 if (qctable) {
166 *qctable = NULL;
167 }
168 /*Check that ilist_fringe and collapse_params are not pointer to NULL
169 * and that there is at least one image in the hdrl imagelist */
170 cpl_error_ensure(ilist_fringe && collapse_params , CPL_ERROR_NULL_INPUT,
171 goto cleanup, "NULL input imagelist or parameter");
172 cpl_error_ensure(hdrl_imagelist_get_size(ilist_fringe) > 0 ,
173 CPL_ERROR_NULL_INPUT, goto cleanup,
174 "input imagelist is empty");
175
176 cpl_size nx, ny, nx1, ny1;
177 nx = hdrl_image_get_size_x(hdrl_imagelist_get_const(ilist_fringe, 0));
178 ny = hdrl_image_get_size_y(hdrl_imagelist_get_const(ilist_fringe, 0));
179
180 /* If there is a ilist_obj check the size and the dimensions: */
181 if (ilist_obj != NULL) {
182 cpl_error_ensure(hdrl_imagelist_get_size(ilist_fringe) ==
183 cpl_imagelist_get_size(ilist_obj),
184 CPL_ERROR_INCOMPATIBLE_INPUT, goto cleanup,
185 "size of fringe and object image list does "
186 "not match");
187
188 nx1 = cpl_image_get_size_x(cpl_imagelist_get_const(ilist_obj, 0));
189 ny1 = cpl_image_get_size_y(cpl_imagelist_get_const(ilist_obj, 0));
190 cpl_error_ensure(nx == nx1, CPL_ERROR_INCOMPATIBLE_INPUT, goto cleanup,
191 "size of fringe image and object mask does not match");
192 cpl_error_ensure(ny == ny1, CPL_ERROR_INCOMPATIBLE_INPUT, goto cleanup,
193 "size of fringe image and object mask does not match");
194 }
195
196 /* If there is a static mask check the dimensions: */
197 if (stat_mask != NULL) {
198 /* If there is a stat_mask check the dimensions: */
199 cpl_error_ensure(nx == cpl_mask_get_size_x(stat_mask),
200 CPL_ERROR_INCOMPATIBLE_INPUT, goto cleanup,
201 "size of fringe image and fringe mask does not match");
202 cpl_error_ensure(ny == cpl_mask_get_size_y(stat_mask),
203 CPL_ERROR_INCOMPATIBLE_INPUT, goto cleanup,
204 "size of fringe image and fringe mask does not match");
205 }
206
207/* This algorithm combines bpm (ilist_fringe), object list (from ilist_obj)
208 * and static mask (stat_mask) for the fringe computation but uses only
209 * the combined bpm and object list for the collapsing */
210
211/* !! This algorithm directly works on the ilist_fringe object thus it
212 * modifies ilist_fringe !! */
213
214 double bkg_level = 0.;
215 double fringe_level = 0.;
216
217 cpl_size isize = hdrl_imagelist_get_size(ilist_fringe);
218 cpl_msg_debug(cpl_func, "Measure fringe amplitudes");
219 if (qctable != NULL) {
220 *qctable = cpl_table_new(isize);
221 cpl_table_new_column(*qctable, "Background_level", CPL_TYPE_DOUBLE);
222 cpl_table_new_column(*qctable, "Fringe_amplitude", CPL_TYPE_DOUBLE);
223 }
224 for (cpl_size i = 0; i < isize; i++) {
225
226 hdrl_image * this_himg = hdrl_imagelist_get(ilist_fringe, i);
227
228 cpl_mask *this_fmsk = cpl_mask_duplicate(hdrl_image_get_mask(this_himg));
229 if (ilist_obj != NULL) {
230 const cpl_image * obj = cpl_imagelist_get_const(ilist_obj, i);
231 cpl_mask * obj_mask = cpl_mask_threshold_image_create(obj,
232 -0.5, 0.5) ;
233 cpl_mask_not(obj_mask) ;
234 cpl_mask_or(this_fmsk, obj_mask);
235 cpl_mask_delete(obj_mask) ;
236 }
237
238 /* Add the object mask to the bad pixel mask for the collapsing */
239 hdrl_image_reject_from_mask(this_himg, this_fmsk);
240
241 /* Add the static mask to the bad pixel mask and object mask for the
242 * amplitude calculation */
243
244 if (stat_mask != NULL) {
245 cpl_mask_or(this_fmsk, stat_mask);
246 }
247 cpl_errorstate prestate = cpl_errorstate_get();
248 cpl_matrix *cur_amplitudes = hdrl_mime_fringe_amplitudes(
249 hdrl_image_get_image_const(this_himg), this_fmsk);
250
251 /* Handle situation where fit is not converging */
252 if (!cpl_errorstate_is_equal(prestate)) {
253 cpl_msg_warning(cpl_func, "Background level and fringe amplitude "
254 "could not be determined! Assuming a background"
255 " level of 0 and a fringe amplitude of 1");
256 bkg_level = 0.;
257 fringe_level = 1.;
258 /* Reset error code */
259 cpl_errorstate_set(prestate);
260 }
261 else {
262 bkg_level = cpl_matrix_get(cur_amplitudes, 0, 0);
263 fringe_level = cpl_matrix_get(cur_amplitudes, 1, 0);
264 }
265
266
267 double fringe_amplitude = fringe_level - bkg_level;
268 if (qctable != NULL) {
269 cpl_table_set_double(*qctable,"Background_level", i, bkg_level);
270 cpl_table_set_double(*qctable,"Fringe_amplitude", i, fringe_amplitude);
271 }
272 cpl_msg_info(cpl_func, "img: %04d Bkg: %12.6g Amplitude: %12.6g",
273 (int)i+1, bkg_level, fringe_amplitude);
274
275 cpl_msg_debug(cpl_func, "Rescaling image");
276
277 hdrl_image_sub_scalar(this_himg, (hdrl_value){bkg_level, 0.});
278
279 hdrl_image_div_scalar(this_himg, (hdrl_value){fringe_amplitude, 0.});
280
281 cpl_matrix_delete(cur_amplitudes);
282 cpl_mask_delete(this_fmsk);
283 }
284
285 cpl_msg_debug(cpl_func, "Combining the normalized fringes generating"
286 " the master-fringe");
287 hdrl_imagelist_collapse(ilist_fringe, collapse_params, master,
288 contrib_map);
289
290cleanup:
291 if (cpl_error_get_code() != CPL_ERROR_NONE) {
292 if (qctable) {
293 cpl_table_delete(*qctable);
294 *qctable = NULL;
295 }
296 if (master) {
297 *master = NULL;
298 }
299 if (contrib_map) {
300 *contrib_map = NULL;
301 }
302 }
303
304 return cpl_error_get_code();
305
306}
307
308/*---------------------------------------------------------------------------*/
344/*---------------------------------------------------------------------------*/
345cpl_error_code
346hdrl_fringe_correct(hdrl_imagelist * ilist_fringe, const cpl_imagelist * ilist_obj,
347 const cpl_mask * stat_mask, const hdrl_image * masterfringe,
348 cpl_table ** qctable)
349{
350 if (qctable) {
351 *qctable = NULL;
352 }
353 /*Check that ilist_fringe and masterfringe are not NULL pointers
354 * and that there is at least one image in the hdrl imagelist */
355 cpl_ensure_code(ilist_fringe && masterfringe, CPL_ERROR_NULL_INPUT);
356 cpl_ensure_code(hdrl_imagelist_get_size(ilist_fringe) > 0 ,
357 CPL_ERROR_NULL_INPUT);
358
359 cpl_size nx, ny, nx1, ny1;
360 nx = hdrl_image_get_size_x(hdrl_imagelist_get_const(ilist_fringe, 0));
361 ny = hdrl_image_get_size_y(hdrl_imagelist_get_const(ilist_fringe, 0));
362
363 /*Check the dimension of the masterfringe image*/
364 nx1 = hdrl_image_get_size_x(masterfringe);
365 ny1 = hdrl_image_get_size_y(masterfringe);
366 cpl_ensure_code(nx == nx1, CPL_ERROR_INCOMPATIBLE_INPUT );
367 cpl_ensure_code(ny == ny1, CPL_ERROR_INCOMPATIBLE_INPUT );
368
369
370 /* If there is a ilist_obj check the size and the dimensions: */
371 if (ilist_obj != NULL) {
372 cpl_ensure_code(hdrl_imagelist_get_size(ilist_fringe) ==
373 cpl_imagelist_get_size(ilist_obj),
374 CPL_ERROR_INCOMPATIBLE_INPUT);
375
376 nx1 = cpl_image_get_size_x(cpl_imagelist_get_const(ilist_obj, 0));
377 ny1 = cpl_image_get_size_y(cpl_imagelist_get_const(ilist_obj, 0));
378 cpl_ensure_code(nx == nx1, CPL_ERROR_INCOMPATIBLE_INPUT );
379 cpl_ensure_code(ny == ny1, CPL_ERROR_INCOMPATIBLE_INPUT );
380 }
381
382 /* If there is a static mask check the dimensions: */
383 if (stat_mask != NULL) {
384 /* If there is a stat_mask check the dimensions: */
385 cpl_ensure_code(nx == cpl_mask_get_size_x(stat_mask),
386 CPL_ERROR_INCOMPATIBLE_INPUT );
387 cpl_ensure_code(ny == cpl_mask_get_size_y(stat_mask),
388 CPL_ERROR_INCOMPATIBLE_INPUT );
389 }
390 /* !! This algorithm directly works on the ilist_fringe object thus it
391 * modifies ilist_fringe !! */
392
393 /* Do the actual fringemap correction */
394 double bkg_level = 0.;
395 double fringe_level = 0.;
396
397 cpl_size isize = hdrl_imagelist_get_size(ilist_fringe);
398 cpl_msg_debug(cpl_func, "Measure fringe amplitudes");
399
400 if (qctable != NULL) {
401 *qctable = cpl_table_new(isize);
402 cpl_table_new_column(*qctable, "Background_level", CPL_TYPE_DOUBLE);
403 cpl_table_new_column(*qctable, "Fringe_amplitude", CPL_TYPE_DOUBLE);
404 }
405 for (cpl_size i = 0; i < isize; i++) {
406
407 hdrl_image * this_himg = hdrl_imagelist_get(ilist_fringe, i);
408 hdrl_image * this_masterfringe = hdrl_image_duplicate(masterfringe);
409
410 cpl_mask *this_fmsk = cpl_mask_duplicate(hdrl_image_get_mask(this_himg));
411 if (stat_mask != NULL) {
412 cpl_mask_or(this_fmsk, stat_mask);
413 }
414
415 if (ilist_obj != NULL) {
416 const cpl_image * obj = cpl_imagelist_get_const(ilist_obj, i);
417 cpl_mask * obj_mask = cpl_mask_threshold_image_create(obj,
418 -0.5, 0.5) ;
419 cpl_mask_not(obj_mask) ;
420 cpl_mask_or(this_fmsk, obj_mask);
421 cpl_mask_delete(obj_mask) ;
422 }
423
424 cpl_errorstate prestate = cpl_errorstate_get();
425
426 cpl_matrix * cur_amplitudes = hdrl_mime_fringe_amplitudes_ls(
427 hdrl_image_get_image_const(this_himg), this_fmsk,
428 hdrl_image_get_image_const(this_masterfringe));
429
430 /* Handle situation where fit is not converging */
431 if (!cpl_errorstate_is_equal(prestate)) {
432 cpl_msg_warning(cpl_func, "Background level and fringe amplitude "
433 "could not be determined! Assuming a background"
434 " level of 0 and a fringe amplitude of 0, i.e. "
435 "no correction will be applied to this image");
436 bkg_level = 0.;
437 fringe_level = 0.;
438 /* Reset error code */
439 cpl_errorstate_set(prestate);
440 }
441 else {
442 bkg_level = cpl_matrix_get(cur_amplitudes, 0, 0);
443 fringe_level = cpl_matrix_get(cur_amplitudes, 1, 0);
444 }
445
446 double fringe_amplitude = fringe_level - bkg_level;
447
448 if (qctable != NULL) {
449 cpl_table_set_double(*qctable,"Background_level", i, bkg_level);
450 cpl_table_set_double(*qctable,"Fringe_amplitude", i, fringe_amplitude);
451 }
452 cpl_msg_info(cpl_func, "img: %04d Bkg: %12.6g Amplitude: %12.6g",
453 (int)i+1, bkg_level, fringe_amplitude);
454
455 cpl_msg_debug(cpl_func, "Rescaling masterfringe");
456 hdrl_image_mul_scalar(this_masterfringe, (hdrl_value){fringe_amplitude, 0.});
457
458 cpl_msg_debug(cpl_func, "Subtract rescaled masterfringe");
459 hdrl_image_sub_image(this_himg, this_masterfringe);
460
461 hdrl_image_delete(this_masterfringe);
462 cpl_matrix_delete(cur_amplitudes);
463 cpl_mask_delete(this_fmsk);
464 }
465
466 if (cpl_error_get_code() != CPL_ERROR_NONE && qctable != NULL) {
467 cpl_table_delete(*qctable);
468 *qctable = NULL;
469 }
470
471 return cpl_error_get_code();
472
473}
474
476/*---------------------------------------------------------------------------*/
498/*---------------------------------------------------------------------------*/
499cpl_matrix *hdrl_mime_fringe_amplitudes(const cpl_image * img0,
500 const cpl_mask * mask0)
501{
502 cpl_matrix *mdata;
503 cpl_matrix *coeffs;
504 cpl_matrix *x;
505 cpl_matrix *hseries;
506 cpl_matrix *amplitudes;
507 cpl_vector *vals;
508 cpl_vector *params;
509 const double *img_data;
510 const cpl_binary *mask_data;
511 double *md;
512 double *par;
513 double mean, stdev, bkg_amp, fringe_amp;
514 int nx, ny, n, size, ns;
515 int msize, i;
516
517 cpl_ensure(img0 != NULL, CPL_ERROR_NULL_INPUT, NULL);
518 cpl_ensure(mask0 != NULL, CPL_ERROR_NULL_INPUT, NULL);
519 cpl_ensure(cpl_image_get_type(img0) == CPL_TYPE_DOUBLE,
520 CPL_ERROR_INVALID_TYPE, NULL);
521
522/* setting parameters
523 n number of the Hermite functions
524 */
525 n = 20;
526
527/* getting image parameters and stats
528 nx, ny dimensions of images
529 */
530 nx = cpl_image_get_size_x(img0);
531 ny = cpl_image_get_size_y(img0);
532 size = nx * ny;
533 msize = size - cpl_mask_count(mask0);
534 /* check that at least some region have been flagged */
535 cpl_ensure(msize > 0 , CPL_ERROR_ILLEGAL_INPUT,NULL);
536
537/* creating masked image */
538 mdata = cpl_matrix_new(msize, 1);
539 md = cpl_matrix_get_data(mdata);
540
541 img_data = cpl_image_get_data_double_const(img0);
542 mask_data = cpl_mask_get_data_const(mask0);
543 for (i = 0; i < size; i++, mask_data++, img_data++)
544 {
545 if (*mask_data == CPL_BINARY_0)
546 {
547 *md = (double) *img_data;
548 md++;
549 }
550 }
551 mean = cpl_matrix_get_mean(mdata);
552 stdev = cpl_matrix_get_stdev(mdata);
553
554/* computing the Hermite coefficients */
555 coeffs = hdrl_mime_hermite_functions_sums_create(n, mean, stdev, mdata);
556 cpl_matrix_multiply_scalar(coeffs, 1.0 / msize);
557
558/* reconstructing the density function as the Hermite series */
559 ns = 1000;
560 x = hdrl_mime_matrix_linspace_create(ns, mean - 4.0 * stdev,
561 mean + 4.0 * stdev);
562 hseries = hdrl_mime_hermite_series_create(n, mean, stdev, coeffs, x);
563
564/* computing the parameters of the Gaussian mixture */
565
566 params = cpl_vector_new(6);
567 par = cpl_vector_get_data(params);
568
569 par[0] = 0.62 / (sqrt(CPL_MATH_PI) * stdev);
570 par[1] = mean - 0.4 * stdev;
571 par[2] = 0.58 * stdev;
572
573 par[3] = 0.57 / (sqrt(CPL_MATH_PI) * stdev);
574 par[4] = mean + 0.3 * stdev;
575 par[5] = 0.61 * stdev;
576
577 vals = cpl_vector_wrap(ns, cpl_matrix_get_data(hseries));
578 cpl_fit_lvmq(x, NULL, vals, NULL, params, NULL, hdrl_mime_gmix1,
579 hdrl_mime_gmix_derivs1, CPL_FIT_LVMQ_TOLERANCE, CPL_FIT_LVMQ_COUNT,
580 CPL_FIT_LVMQ_MAXITER, NULL, NULL, NULL);
581
582 bkg_amp = (par[1] > par[4] ? par[4] : par[1]);
583 fringe_amp = (par[1] > par[4] ? par[1] : par[4]);
584 amplitudes = cpl_matrix_new(2, 1);
585 cpl_matrix_set(amplitudes, 0, 0, bkg_amp);
586 cpl_matrix_set(amplitudes, 1, 0, fringe_amp);
587
588/* cleaning up */
589 cpl_matrix_delete(mdata);
590 cpl_matrix_delete(coeffs);
591 cpl_matrix_delete(x);
592 cpl_matrix_delete(hseries);
593 cpl_vector_unwrap(vals);
594 cpl_vector_delete(params);
595
596 return amplitudes;
597}
598
599
600/*---------------------------------------------------------------------------*/
619/*---------------------------------------------------------------------------*/
620cpl_matrix *hdrl_mime_fringe_amplitudes_ls(const cpl_image * img0,
621 const cpl_mask * mask0, const cpl_image * fringe0)
622{
623 cpl_matrix *mdata;
624 cpl_matrix *fdata;
625 cpl_matrix *cols12;
626 cpl_matrix *coeffs;
627 cpl_matrix *amplitudes;
628 const double *img_data;
629 const cpl_binary *mask_data;
630 const double *fringe_data;
631 double *md;
632 double *fd;
633 int nx, ny, size;
634 int msize, i;
635
636 cpl_ensure(img0 != NULL, CPL_ERROR_NULL_INPUT, NULL);
637 cpl_ensure(mask0 != NULL, CPL_ERROR_NULL_INPUT, NULL);
638 cpl_ensure(fringe0 != NULL, CPL_ERROR_NULL_INPUT, NULL);
639 cpl_ensure(cpl_image_get_type(img0) == CPL_TYPE_DOUBLE,
640 CPL_ERROR_INVALID_TYPE, NULL);
641 cpl_ensure(cpl_image_get_type(fringe0) == CPL_TYPE_DOUBLE,
642 CPL_ERROR_INVALID_TYPE, NULL);
643
644/* getting image parameters and stats
645 nx, ny dimensions of images
646 */
647 nx = cpl_image_get_size_x(img0);
648 ny = cpl_image_get_size_y(img0);
649 size = nx * ny;
650 msize = size - cpl_mask_count(mask0);
651 /* check that at least some region have been flagged */
652 cpl_ensure(msize > 0 , CPL_ERROR_ILLEGAL_INPUT,NULL);
653
654/* creating masked image and masked fringe */
655 mdata = cpl_matrix_new(msize, 1);
656 md = cpl_matrix_get_data(mdata);
657
658 fdata = cpl_matrix_new(msize, 1);
659 fd = cpl_matrix_get_data(fdata);
660
661 img_data = cpl_image_get_data_double_const(img0);
662 mask_data = cpl_mask_get_data_const(mask0);
663 fringe_data = cpl_image_get_data_double_const(fringe0);
664 for (i = 0; i < size; i++, img_data++, mask_data++, fringe_data++)
665 {
666 if (*mask_data == CPL_BINARY_0)
667 {
668 *md = (double) *img_data;
669 *fd = (double) *fringe_data;
670 md++;
671 fd++;
672 }
673 }
674
675/* computing the least squares coefficients */
676 cols12 = cpl_matrix_new(msize, 2);
677 cpl_matrix_fill(cols12, 1.0);
678 cpl_matrix_copy(cols12, fdata, 0, 0);
679 coeffs = hdrl_mime_linalg_solve_tikhonov(cols12, mdata, 1.0e-10);
680
681/* computing the amplitudes */
682 amplitudes = cpl_matrix_new(2, 1);
683 cpl_matrix_set(amplitudes, 0, 0, cpl_matrix_get(coeffs, 1, 0));
684 cpl_matrix_set(amplitudes, 1, 0, cpl_matrix_get(coeffs, 0, 0) +
685 cpl_matrix_get(coeffs, 1, 0));
686
687/* cleaning up */
688 cpl_matrix_delete(mdata);
689 cpl_matrix_delete(fdata);
690 cpl_matrix_delete(cols12);
691 cpl_matrix_delete(coeffs);
692
693 return amplitudes;
694}
695
696
697/*---------------------------------------------------------------------------*/
708/*---------------------------------------------------------------------------*/
709int hdrl_mime_gmix_derivs1(const double x[], const double params[],
710 double result[])
711{
712 double a1, m1, sigma1;
713 double a2, m2, sigma2;
714 double tmp;
715
716/* initializing */
717 a1 = params[0];
718 m1 = params[1];
719 sigma1 = params[2];
720
721 a2 = params[3];
722 m2 = params[4];
723 sigma2 = params[5];
724
725/* evaluating */
726 tmp = (x[0] - m1) / sigma1;
727 result[0] = exp(-0.5 * tmp * tmp);
728 result[1] = a1 * exp(-0.5 * tmp * tmp);
729 result[1] *= tmp / sigma1;
730 result[2] = a1 * exp(-0.5 * tmp * tmp);
731 result[2] *= tmp * tmp / sigma1;
732
733 tmp = (x[0] - m2) / sigma2;
734 result[3] = exp(-0.5 * tmp * tmp);
735 result[4] = a2 * exp(-0.5 * tmp * tmp);
736 result[4] *= tmp / sigma2;
737 result[5] = a2 * exp(-0.5 * tmp * tmp);
738 result[5] *= tmp * tmp / sigma2;
739
740 return 0;
741}
742
743/*---------------------------------------------------------------------------*/
754/*---------------------------------------------------------------------------*/
755int hdrl_mime_gmix1(const double x[], const double params[], double *result)
756{
757 double a1, m1, sigma1;
758 double a2, m2, sigma2;
759 double tmp;
760
761/* initializing */
762 a1 = params[0];
763 m1 = params[1];
764 sigma1 = params[2];
765
766 a2 = params[3];
767 m2 = params[4];
768 sigma2 = params[5];
769
770/* evaluating */
771 tmp = (x[0] - m1) / sigma1;
772 result[0] = a1 * exp(-0.5 * tmp * tmp);
773 tmp = (x[0] - m2) / sigma2;
774 result[0] += a2 * exp(-0.5 * tmp * tmp);
775
776 return 0;
777}
778
779/*---------------------------------------------------------------------------*/
795/*---------------------------------------------------------------------------*/
796cpl_matrix *hdrl_mime_hermite_series_create(int n, double center,
797 double scale, const cpl_matrix * coeffs, const cpl_matrix * x)
798{
799 cpl_matrix *series;
800 double *ms;
801 const double *mc;
802 const double *mx;
803 double rt;
804 int i, k, size;
805
806/* testing input */
807 if (x == NULL || coeffs == NULL)
808 {
809 cpl_error_set(cpl_func, CPL_ERROR_NULL_INPUT);
810 return NULL;
811 }
812
813 if (n < 1 || scale <= 0.0)
814 {
815 cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_INPUT);
816 return NULL;
817 }
818
819/* The specific dimensions of the matrix x are not used, only its size. */
820 size = cpl_matrix_get_nrow(x) * cpl_matrix_get_ncol(x);
821 mx = cpl_matrix_get_data_const(x);
822 mc = cpl_matrix_get_data_const(coeffs);
823
824/* allocating memory */
825 series = cpl_matrix_new(size, 1);
826 ms = cpl_matrix_get_data(series);
827
828/* computing the normalization constant */
829 rt = 1.0 / sqrt(sqrt(CPL_MATH_PI));
830
831/* evaluating the Hermite functions at the samples */
832 for (k = 0; k < size; k++)
833 {
834 double xk = (mx[k] - center) / scale;
835 double tmp1 = rt * exp(-0.5 * xk * xk);
836 double tmp2 = rt * sqrt(2.0) * xk * exp(-0.5 * xk * xk);
837 for (i = 2; i < n + 2; i++)
838 {
839 double tmp3 = sqrt(2) * xk * tmp2 - sqrt(i - 1) * tmp1;
840 tmp3 = tmp3 / sqrt(i);
841
842 ms[k] += tmp1 * mc[i - 2];
843 tmp1 = tmp2;
844 tmp2 = tmp3;
845 }
846 }
847
848/* normalizing */
849 cpl_matrix_multiply_scalar(series, 1 / sqrt(scale));
850
851 return series;
852}
853
854
855
856/*---------------------------------------------------------------------------*/
872/*---------------------------------------------------------------------------*/
873cpl_matrix *hdrl_mime_hermite_functions_sums_create(int n, double center,
874 double scale, const cpl_matrix * x)
875{
876 cpl_matrix *sums;
877 double *ms;
878 const double *mx;
879 double rt;
880 int i, k, size;
881 double tmp3;
882 double sqrt_i[n + 2];
883 double rsqrt_i[n + 2];
884
885/* testing input */
886 if (x == NULL)
887 {
888 cpl_error_set(cpl_func, CPL_ERROR_NULL_INPUT);
889 return NULL;
890 }
891
892 if (n < 1 || scale <= 0.0)
893 {
894 cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_INPUT);
895 return NULL;
896 }
897
898/* The specific dimensions of the matrix x are not used, only its size. */
899 size = cpl_matrix_get_nrow(x) * cpl_matrix_get_ncol(x);
900 mx = cpl_matrix_get_data_const(x);
901
902/* allocating memory */
903 sums = cpl_matrix_new(n, 1);
904 ms = cpl_matrix_get_data(sums);
905
906/* computing the normalization constant */
907 rt = 1.0 / sqrt(sqrt(CPL_MATH_PI));
908 for (i = 1; i < n + 2; i++) {
909 sqrt_i[i] = sqrt(i);
910 rsqrt_i[i] = 1. / sqrt_i[i];
911 }
912
913/* evaluating the Hermite functions at the samples */
914 for (k = 0; k < size; k++)
915 {
916 double xk = (mx[k] - center) / scale;
917 double tmp1 = rt * exp(-0.5 * xk * xk);
918 double tmp2 = rt * sqrt(2.0) * xk * exp(-0.5 * xk * xk);
919 for (i = 2; i < n + 2; i++)
920 {
921 tmp3 = sqrt(2) * xk * tmp2 - sqrt_i[i - 1] * tmp1;
922 tmp3 = tmp3 * rsqrt_i[i];
923
924 ms[i - 2] += tmp1;
925 tmp1 = tmp2;
926 tmp2 = tmp3;
927 }
928 }
929
930/* normalizing */
931 cpl_matrix_multiply_scalar(sums, 1 / sqrt(scale));
932
933 return sums;
934}
cpl_error_code hdrl_fringe_correct(hdrl_imagelist *ilist_fringe, const cpl_imagelist *ilist_obj, const cpl_mask *stat_mask, const hdrl_image *masterfringe, cpl_table **qctable)
Scales and subtracts the master fringe from the images.
Definition: hdrl_fringe.c:346
cpl_error_code hdrl_fringe_compute(hdrl_imagelist *ilist_fringe, const cpl_imagelist *ilist_obj, const cpl_mask *stat_mask, const hdrl_parameter *collapse_params, hdrl_image **master, cpl_image **contrib_map, cpl_table **qctable)
Calculates the master fringe and contribution map based on the.
Definition: hdrl_fringe.c:160
cpl_error_code hdrl_image_sub_image(hdrl_image *self, const hdrl_image *other)
Subtract two images, store the result in the first image.
cpl_error_code hdrl_image_reject_from_mask(hdrl_image *self, const cpl_mask *map)
set bpm of hdrl_image
Definition: hdrl_image.c:407
cpl_error_code hdrl_image_div_scalar(hdrl_image *self, hdrl_value value)
Elementwise division of an image with a scalar.
cpl_error_code hdrl_image_mul_scalar(hdrl_image *self, hdrl_value value)
Elementwise multiplication of an image with a scalar.
hdrl_image * hdrl_image_duplicate(const hdrl_image *himg)
copy hdrl_image
Definition: hdrl_image.c:391
cpl_mask * hdrl_image_get_mask(hdrl_image *himg)
get cpl bad pixel mask from image
Definition: hdrl_image.c:157
cpl_size hdrl_image_get_size_y(const hdrl_image *self)
return size of Y dimension of image
Definition: hdrl_image.c:540
cpl_size hdrl_image_get_size_x(const hdrl_image *self)
return size of X dimension of image
Definition: hdrl_image.c:525
const cpl_image * hdrl_image_get_image_const(const hdrl_image *himg)
get data as cpl image
Definition: hdrl_image.c:118
void hdrl_image_delete(hdrl_image *himg)
delete hdrl_image
Definition: hdrl_image.c:379
cpl_error_code hdrl_image_sub_scalar(hdrl_image *self, hdrl_value value)
Elementwise subtraction of a scalar from an image.
const hdrl_image * hdrl_imagelist_get_const(const hdrl_imagelist *himlist, cpl_size inum)
Get an image from a list of images.
cpl_size hdrl_imagelist_get_size(const hdrl_imagelist *himlist)
Get the number of images in the imagelist.
cpl_error_code hdrl_imagelist_collapse(const hdrl_imagelist *himlist, const hdrl_parameter *param, hdrl_image **out, cpl_image **contrib)
collapsing of image list
hdrl_image * hdrl_imagelist_get(const hdrl_imagelist *himlist, cpl_size inum)
Get an image from a list of images.
cpl_matrix * hdrl_mime_linalg_solve_tikhonov(const cpl_matrix *mat, const cpl_matrix *rhs, double alpha)
Solve an overdetermined linear system in the least-squares sense.
cpl_matrix * hdrl_mime_matrix_linspace_create(int n, double a, double b)
Create equally spaced nodes.