29#define REGDEBUG_FULL 0
35#include "hdrl_lacosmics.h"
36#include "hdrl_utils.h"
85} hdrl_lacosmic_parameter;
88static hdrl_parameter_typeobj hdrl_lacosmic_parameter_type = {
89 HDRL_PARAMETER_LACOSMIC,
90 (hdrl_alloc *)&cpl_malloc,
91 (hdrl_free *)&cpl_free,
93 sizeof(hdrl_lacosmic_parameter),
117 hdrl_lacosmic_parameter * p = (hdrl_lacosmic_parameter *)
118 hdrl_parameter_new(&hdrl_lacosmic_parameter_type);
119 p->sigma_lim = sigma_lim ;
121 p->max_iter = max_iter ;
122 return (hdrl_parameter *)p;
133 const hdrl_parameter * param)
135 const hdrl_lacosmic_parameter * param_loc = (
const hdrl_lacosmic_parameter *)param;
137 cpl_error_ensure(param != NULL, CPL_ERROR_NULL_INPUT,
138 return CPL_ERROR_NULL_INPUT,
"NULL Input Parameters");
140 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
141 "Expected LaCosmic parameter") ;
143 cpl_error_ensure(param_loc->max_iter > 0, CPL_ERROR_ILLEGAL_INPUT,
144 return CPL_ERROR_ILLEGAL_INPUT,
"max_iter must be >0");
145 cpl_error_ensure(param_loc->f_lim >= 0, CPL_ERROR_ILLEGAL_INPUT,
146 return CPL_ERROR_ILLEGAL_INPUT,
"f_lim must be >=0");
147 cpl_error_ensure(param_loc->sigma_lim >= 0, CPL_ERROR_ILLEGAL_INPUT,
148 return CPL_ERROR_ILLEGAL_INPUT,
"sigma_lim must be >=0");
150 return CPL_ERROR_NONE ;
162 return hdrl_parameter_check_type(self, &hdrl_lacosmic_parameter_type);
173 const hdrl_parameter * p)
175 cpl_ensure(p, CPL_ERROR_NULL_INPUT, -1.0);
176 return ((
const hdrl_lacosmic_parameter *)p)->sigma_lim;
187 const hdrl_parameter * p)
189 cpl_ensure(p, CPL_ERROR_NULL_INPUT, -1.0);
190 return ((
const hdrl_lacosmic_parameter *)p)->f_lim;
201 const hdrl_parameter * p)
203 cpl_ensure(p, CPL_ERROR_NULL_INPUT, -1);
204 return ((
const hdrl_lacosmic_parameter *)p)->max_iter;
223 const char *base_context,
225 const hdrl_parameter *defaults)
227 cpl_ensure(prefix && base_context && defaults,
228 CPL_ERROR_NULL_INPUT, NULL);
231 CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
233 cpl_parameterlist *parlist = cpl_parameterlist_new();
240 hdrl_setup_vparameter(parlist, prefix,
".",
"",
"sigma_lim", base_context,
241 "Poisson fluctuation threshold to flag cosmics"
242 "(see van Dokkum, PASP,113,2001,p1420-27).",
243 CPL_TYPE_DOUBLE, sigma_lim_def) ;
246 hdrl_setup_vparameter(parlist, prefix,
".",
"",
"f_lim", base_context,
247 "Minimum contrast between the Laplacian image and the fine "
248 "structure image that a point must have to be flagged as cosmics"
249 , CPL_TYPE_DOUBLE, f_lim_def) ;
252 hdrl_setup_vparameter(parlist, prefix,
".",
"",
"max_iter", base_context,
253 "Maximum number of alghoritm iterations", CPL_TYPE_INT, max_iter_def) ;
255 if (cpl_error_get_code()) {
256 cpl_parameterlist_delete(parlist);
277 const cpl_parameterlist * parlist,
280 cpl_ensure(prefix && parlist, CPL_ERROR_NULL_INPUT, NULL);
282 const cpl_parameter * par;
283 double sigma_lim, f_lim;
288 par=cpl_parameterlist_find_const(parlist, name);
289 sigma_lim = cpl_parameter_get_double(par);
294 par=cpl_parameterlist_find_const(parlist, name);
295 f_lim = cpl_parameter_get_double(par);
300 par=cpl_parameterlist_find_const(parlist, name);
301 max_iter = cpl_parameter_get_int(par);
304 if (cpl_error_get_code()) {
305 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
306 "Error while parsing parameterlist with prefix %s", prefix);
349 const hdrl_image * ima_in,
350 const hdrl_parameter * params)
352 cpl_image * sci_data;
353 cpl_image * sci_error = NULL;
354 cpl_mask * sci_mask = NULL;
355 cpl_image * laplacian_redu_data = NULL;
356 cpl_image * subs2_data = NULL;
357 cpl_image * s_data = NULL;
358 cpl_image * f_data = NULL;
359 cpl_image * r_data = NULL;
362 intptr_t subs2_nx = 0;
363 intptr_t subs2_ny = 0;
365 cpl_binary * psci_mask = NULL;
366 cpl_binary * pout_mask = NULL;
368 cpl_matrix * laplacian_kernel = NULL;
370 cpl_mask * median3_kernel = NULL;
371 cpl_mask * median5_kernel = NULL;
372 cpl_mask * median7_kernel = NULL;
373 cpl_mask * lastiter_mask = NULL;
374 cpl_mask * out_mask = NULL;
376 double * psci_data = NULL;
377 double * psci_error = NULL;
378 double * psubs2_data = NULL;
379 double * plaplacian_data = NULL;
380 double * plaplacian_redu_data = NULL;
381 double * psci_median3_data = NULL;
382 double * ps_data = NULL;
383 double * pf_data = NULL;
384 double * pr_data = NULL;
392 cpl_ensure(ima_in, CPL_ERROR_NULL_INPUT, NULL);
394 CPL_ERROR_ILLEGAL_INPUT, NULL) ;
397 CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
399 CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
402 const hdrl_lacosmic_parameter * p_loc = (
const hdrl_lacosmic_parameter *)params ;
409 sci_mask = cpl_mask_new(cpl_image_get_size_x(sci_data),
410 cpl_image_get_size_y(sci_data));
416 laplacian_kernel = cpl_matrix_new(3, 3);
417 cpl_matrix_set(laplacian_kernel, 0, 0, 0.0);
418 cpl_matrix_set(laplacian_kernel, 0, 1, -1.0);
419 cpl_matrix_set(laplacian_kernel, 0, 2, 0.0);
420 cpl_matrix_set(laplacian_kernel, 1, 0, -1.0);
421 cpl_matrix_set(laplacian_kernel, 1, 1, 4.0);
422 cpl_matrix_set(laplacian_kernel, 1, 2, -1.0);
423 cpl_matrix_set(laplacian_kernel, 2, 0, 0.0);
424 cpl_matrix_set(laplacian_kernel, 2, 1, -1.0);
425 cpl_matrix_set(laplacian_kernel, 2, 2, 0.0);
428 median3_kernel = cpl_mask_new(3, 3);
429 cpl_mask_not(median3_kernel);
432 median5_kernel = cpl_mask_new(5, 5);
433 cpl_mask_not(median5_kernel);
436 median7_kernel = cpl_mask_new(7, 7);
437 cpl_mask_not(median7_kernel);
439 out_mask = cpl_mask_new(cpl_mask_get_size_x(sci_mask),
440 cpl_mask_get_size_y(sci_mask));
442 nx = cpl_image_get_size_x(sci_data);
443 ny = cpl_image_get_size_y(sci_data);
446 psci_data = cpl_image_get_data_double(sci_data);
447 psci_error = cpl_image_get_data_double(sci_error);
449 psci_mask = cpl_mask_get_data(sci_mask);
450 pout_mask = cpl_mask_get_data(out_mask);
455 subs2_data = cpl_image_new(subs2_nx, subs2_ny, CPL_TYPE_DOUBLE);
456 psubs2_data = cpl_image_get_data_double( subs2_data);
458 laplacian_redu_data = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
459 plaplacian_redu_data = cpl_image_get_data_double( laplacian_redu_data);
461 s_data = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
462 ps_data = cpl_image_get_data_double( s_data);
464 f_data = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
465 pf_data = cpl_image_get_data_double( f_data);
467 r_data = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
468 pr_data = cpl_image_get_data_double( r_data);
472 lastiter_mask = cpl_mask_duplicate(out_mask);
476 while (nbiter <= p_loc->max_iter) {
478 double * psci_median3_7_data;
479 double * ps_median_data;
483 if (nbiter > 1 && !hdrl_check_maskequality(lastiter_mask, out_mask)) {
484 cpl_msg_debug(cpl_func,
"Detections of iteration %d and %d are "
485 "identical - stopping here", nbiter-1, nbiter );
488 cpl_mask_delete(lastiter_mask);
489 lastiter_mask = cpl_mask_duplicate(out_mask);
498 for (intptr_t j = 0; j < ny; j++) {
499 intptr_t j_2_subs2_nx = j * 2 * subs2_nx;
500 intptr_t j_nx = j * nx;
501 for (intptr_t i = 0; i < nx; i++) {
502 double val = psci_data[i + j_nx];
503 intptr_t pix = i * 2 + j_2_subs2_nx;
504 psubs2_data[pix] = val;
505 psubs2_data[pix + subs2_nx] = val;
507 psubs2_data[pix] = val;
508 psubs2_data[pix + subs2_nx] = val;
518 cpl_image * laplacian_data = hdrl_parallel_filter_image(subs2_data,
519 laplacian_kernel, NULL,
522 plaplacian_data = cpl_image_get_data_double(laplacian_data);
523 for (intptr_t i = 0; i < subs2_nx * subs2_ny; i++) {
524 if (plaplacian_data[i] < 0.0) {
525 plaplacian_data[i] = 0.0;
528 plaplacian_data[i] *= 8.0;
534 for (intptr_t i = 0; i < subs2_ny; i++) {
535 plaplacian_data[i * subs2_nx + 0] = plaplacian_data[i * subs2_nx + 1];
536 plaplacian_data[i * subs2_nx + subs2_nx - 1] =
537 plaplacian_data[i * subs2_nx + + subs2_nx - 2];
539 for (intptr_t j = 0; j < subs2_nx; j++) {
540 plaplacian_data[0 * subs2_nx + j] = plaplacian_data[1 * subs2_nx + j];
541 plaplacian_data[(subs2_ny - 1) * subs2_nx + j] =
542 plaplacian_data[(subs2_ny - 2) * subs2_nx + j];
546 cpl_image_save(laplacian_data,
"Lpositive.fits", CPL_BPP_IEEE_DOUBLE,
547 NULL, CPL_IO_DEFAULT);
556HDRL_OMP(omp parallel
for)
557 for (intptr_t j = 0; j < ny; j++) {
558 intptr_t j_2_subs2_nx = j * 2 * subs2_nx;
559 intptr_t j_nx = j * nx;
560 for (intptr_t i = 0; i < nx; i++) {
561 intptr_t pix = i * 2 + j_2_subs2_nx;
562 plaplacian_redu_data[i + j_nx] =
563 (plaplacian_data[pix] +
564 plaplacian_data[pix + 1] +
565 plaplacian_data[pix + subs2_nx] +
566 plaplacian_data[pix + subs2_nx + 1]) * 0.25;
569 0.5 * plaplacian_redu_data[i + j_nx] / psci_error[i + j_nx];
574 cpl_image_save(laplacian_redu_data,
"Lplus.fits", CPL_BPP_IEEE_DOUBLE,
575 NULL, CPL_IO_DEFAULT);
579 cpl_image * s_median_data = hdrl_parallel_filter_image(s_data, NULL,
582 ps_median_data = cpl_image_get_data_double( s_median_data);
586 for (intptr_t i = 0; i < nx * ny; i++) {
587 ps_data[i] -= ps_median_data[i];
591 cpl_image_save( s_data,
"S2.fits", CPL_BPP_IEEE_DOUBLE, NULL,
596 cpl_image * sci_median3_data =
597 hdrl_parallel_filter_image(sci_data, NULL,
600 psci_median3_data = cpl_image_get_data_double(sci_median3_data);
604 cpl_image * sci_median3_7_data =
605 hdrl_parallel_filter_image(sci_median3_data, NULL,
608 psci_median3_7_data = cpl_image_get_data_double(sci_median3_7_data);
611 for (intptr_t i = 0; i < nx * ny; i++) {
612 pf_data[i] = psci_median3_data[i] - psci_median3_7_data[i];
618 if (pf_data[i] < 0.01) {
624 cpl_image_save( f_data,
"F.fits", CPL_BPP_IEEE_DOUBLE, NULL,
629 for (intptr_t i = 0; i < nx * ny; i++) {
630 pr_data[i] = plaplacian_redu_data[i] / pf_data[i];
634 cpl_image_save( r_data,
"R.fits", CPL_BPP_IEEE_DOUBLE, NULL,
639 median = cpl_vector_new(24);
640 for (intptr_t j = 0; j < ny - 1; j++) {
641 intptr_t j_nx = j * nx;
642 for (intptr_t i = 0; i < nx - 1; i++) {
643 intptr_t i_plus_j_nx = i + j_nx;
644 if (ps_data[i_plus_j_nx] > p_loc->sigma_lim &&
645 pr_data[i_plus_j_nx] > p_loc->f_lim &&
647 psci_mask[i_plus_j_nx] == CPL_BINARY_0) {
649 cpl_vector* med_vect = NULL;
650 intptr_t li, lj, ui, uj, m;
652 pout_mask[i_plus_j_nx] = CPL_BINARY_1;
653 cpl_msg_debug(cpl_func,
"Detection found at x=%zd y=%zd with "
654 "value=%g", i+1, j+1, psci_data[i_plus_j_nx]);
659 li = CX_MAX(i - 2, 0);
660 lj = CX_MAX(j - 2, 0);
661 ui = CX_MIN(nx, i + 3);
662 uj = CX_MIN(ny, j + 3);
664 for (intptr_t k = lj; k < uj; k++) {
665 intptr_t k_nx = k * nx;
666 for (intptr_t l = li; l < ui; l++) {
667 intptr_t l_plus_k_nx = l + k_nx;
668 if (ps_data[l_plus_k_nx] <= p_loc->sigma_lim &&
670 psci_mask[l_plus_k_nx] == CPL_BINARY_0) {
671 cpl_vector_set(median,m,psci_data[l_plus_k_nx]);
675 if (pr_data[l_plus_k_nx] <= p_loc->f_lim &&
677 psci_mask[l_plus_k_nx] == CPL_BINARY_0) {
678 cpl_vector_set(median,m,psci_data[l_plus_k_nx]);
688 data = cpl_vector_get_data( median);
689 med_vect = cpl_vector_wrap( m, data);
690 psci_data[i_plus_j_nx] = cpl_vector_get_median(med_vect);
691 cpl_msg_debug(cpl_func,
"Detection replaced with value=%g",
692 psci_data[i_plus_j_nx]);
693 cpl_vector_unwrap(med_vect);
698 cpl_vector_delete(median); median = NULL;
701 cpl_image_delete(laplacian_data);
702 cpl_image_delete(sci_median3_7_data);
703 cpl_image_delete(sci_median3_data);
704 cpl_image_delete(s_median_data);
708 cpl_mask_save( *out_mask,
"CRH_SINGLE.fits", NULL, CPL_IO_DEFAULT);
712 cpl_matrix_delete(laplacian_kernel);
713 cpl_mask_delete(median3_kernel);
714 cpl_mask_delete(median5_kernel);
715 cpl_mask_delete(median7_kernel);
716 cpl_mask_delete(lastiter_mask);
717 cpl_image_delete(laplacian_redu_data);
718 cpl_image_delete(subs2_data);
719 cpl_image_delete(s_data);
720 cpl_image_delete(f_data);
721 cpl_image_delete(r_data);
722 cpl_image_delete(sci_data);
723 cpl_image_delete(sci_error);
724 cpl_mask_delete(sci_mask);
cpl_size hdrl_image_get_size_y(const hdrl_image *self)
return size of Y dimension of image
const cpl_mask * hdrl_image_get_mask_const(const hdrl_image *himg)
get cpl bad pixel mask from image
cpl_size hdrl_image_get_size_x(const hdrl_image *self)
return size of X dimension of image
const cpl_image * hdrl_image_get_error_const(const hdrl_image *himg)
get error as cpl image
const cpl_image * hdrl_image_get_image_const(const hdrl_image *himg)
get data as cpl image
cpl_parameterlist * hdrl_lacosmic_parameter_create_parlist(const char *base_context, const char *prefix, const hdrl_parameter *defaults)
Create parameter list for the LaCosmic computation.
double hdrl_lacosmic_parameter_get_sigma_lim(const hdrl_parameter *p)
Access the sigma_lim in the LaCosmic parameter.
cpl_mask * hdrl_lacosmic_edgedetect(const hdrl_image *ima_in, const hdrl_parameter *params)
Detect bad-pixels / cosmic-rays on a single image.
int hdrl_lacosmic_parameter_get_max_iter(const hdrl_parameter *p)
Access the max_iter in the LaCosmic parameter.
cpl_error_code hdrl_lacosmic_parameter_verify(const hdrl_parameter *param)
Verify basic correctness of the LaCosmic parameters.
hdrl_parameter * hdrl_lacosmic_parameter_parse_parlist(const cpl_parameterlist *parlist, const char *prefix)
Parse parameterlist to create input parameters for the LaCosmic.
double hdrl_lacosmic_parameter_get_f_lim(const hdrl_parameter *p)
Access the f_lim in the LaCosmic parameter.
cpl_boolean hdrl_lacosmic_parameter_check(const hdrl_parameter *self)
Check that the parameter is an LaCosmic parameter.
hdrl_parameter * hdrl_lacosmic_parameter_create(double sigma_lim, double f_lim, int max_iter)
Creates LaCosmic parameters object.
char * hdrl_join_string(const char *sep_, int n,...)
join strings together