ERIS Pipeline Reference Manual 1.9.3
eris_nix_img_supersky_static.c
1/* $Id$
2 *
3 * This file is part of the ERIS/NIX Pipeline
4 * Copyright (C) 2025 European Southern Observatory
5 *
6 * This program is free software; you can redistribute it and/or modify
7 * it under the terms of the GNU General Public License as published by
8 * the Free Software Foundation; either version 2 of the License, or
9 * (at your option) any later version.
10 */
11
12#ifdef HAVE_CONFIG_H
13#include <config.h>
14#endif
15
16/*-----------------------------------------------------------------------------
17 Includes
18 -----------------------------------------------------------------------------*/
19#include <cpl.h>
20#include <hdrl.h>
21#include <math.h>
22#include <string.h>
23
24#include "eris_nix_img_supersky.h"
25#include "eris_nix_utils.h"
26#include "eris_nix_dfs.h"
27#include "eris_pfits.h"
28#include "eris_utils.h"
29/*-----------------------------------------------------------------------------
30 Private function prototypes
31 -----------------------------------------------------------------------------*/
32static hdrl_image * eris_create_sky_flat_first(
33 hdrl_imagelist * image_list,
34 cpl_mask * masks,
35 const char * method,
36 const double kappa,
37 const int niter,
38 const cpl_size plane_id);
39
40static hdrl_image * eris_create_sky_flat_first2(
41 located_imagelist * jitters,
42 located_imagelist ** sky_1st,
43 cpl_mask * mask,
44 const char * method,
45 const double kappa,
46 const int niter,
47 const cpl_size plane_id);
48
49static hdrl_image * eris_create_sky_flat_final(
50 hdrl_imagelist * image_list,
51 cpl_mask ** masks,
52 const char * method,
53 const double kappa,
54 const int niter,
55 const cpl_size plane_id,
56 cpl_image** contrib);
57
58static hdrl_image * eris_create_sky_flat_final2(
59 located_imagelist * jitters,
60 located_imagelist ** sky_2nd,
61 cpl_mask ** masks,
62 const char * method,
63 const double kappa,
64 const int niter,
65 const cpl_size plane_id);
66
67static cpl_mask * eris_detect_sources_hdrl_catalogue(
68 const cpl_image * image,
69 hdrl_parameter* p);
70
71static cpl_mask * eris_dilate_mask(
72 const cpl_mask * input_mask,
73 int radius);
74
75static cpl_error_code
76eris_image_stats(const hdrl_image *in_image,
77 double l_sig,
78 double u_sig,
79 const cpl_mask *mask,
80 double trim_fraction,
81 double *median,
82 double *mean,
83 double *stdev,
84 double *tmean,
85 double *tstd,
86 double *mad,
87 double *min_val,
88 double *max_val);
89
90static cpl_error_code
91eris_first_sky_sub(const hdrl_imagelist* science_images, const hdrl_image* sky_flat_0,
92 hdrl_imagelist** skysub_images);
93
94static cpl_error_code
95eris_locate_and_mask_sources(hdrl_imagelist* skysub_images, const cpl_image* bpm_image,
96 const cpl_parameterlist* parlist, const char* prefix, cpl_mask** source_masks);
97
98
99static cpl_error_code
100eris_locate_and_mask_sources2( located_imagelist * jitters, const cpl_image* bpm_image,
101 const cpl_parameterlist* parlist, const char* prefix, cpl_mask** source_masks);
102
103static cpl_error_code
104eris_load_data_and_crea_error(const char* fname, const cpl_size kk, const cpl_mask* bad_pixel_map,
105 hdrl_imagelist** science_images);
106
107static cpl_error_code
108eris_load_data_and_error_simple(const char* fname, const cpl_size kk, const cpl_mask* bad_pixel_map,
109 hdrl_imagelist** science_images);
110
111static cpl_error_code
112eris_load_and_error_cube(const char* fname, const cpl_mask* bad_pixel_map, cpl_size* kkk,
113 hdrl_imagelist** science_images);
114
115static cpl_error_code
116eris_load_and_error_cube2(const char* fname, const cpl_mask* bad_pixel_map,
117 const char* method, const double kappa, const int niter,
118 hdrl_imagelist** science_images);
119
120static cpl_error_code
121eris_create_final_skysub_products(const hdrl_imagelist* science_iml, const hdrl_image* sky_flat, cpl_frameset** products);
122
123static hdrl_imagelist*
124eris_crea_imagelist_all(cpl_frameset* raw_frames, cpl_mask* bad_pixel_map);
125
126static cpl_boolean
127eris_check_format_is_cube_of_same_size(const cpl_frameset * raw_frames);
128
129static cpl_error_code
130eris_load_and_error_cube_slice(const char* fname, const cpl_mask* bad_pixel_map, const cpl_size k,
131 const cpl_size kkk, hdrl_imagelist** science_images);
132
133
134static cpl_error_code
135eris_load_and_error_cube_robert(const char* fname, const cpl_mask* bad_pixel_map,
136 const char* method, const double kappa, const int niter, const cpl_size kkk,
137 hdrl_imagelist** science_images);
138
139static cpl_error_code
140eris_save_sky_flat_final( const cpl_frameset* raw_frames, const hdrl_image* sky_flat_final, cpl_frameset** product_frames);
141
142static cpl_error_code
143eris_handle_simple_format_case(cpl_frameset * frameset,
144 const cpl_parameterlist * parlist,
145 cpl_frameset * raw_frames,
146 const char * recipe_name,
147 const char* combine_method,
148 const double sigma_clip,
149 const int max_iter,
150 const cpl_boolean debug_data,
151 cpl_mask * bad_pixel_map,
152 cpl_frameset ** product_frames);
153
154
155static cpl_error_code
156eris_handle_simple_format_case2(cpl_frameset * frameset,
157 const cpl_parameterlist * parlist,
158 const char * recipe_name,
159 const char* combine_method,
160 const double sigma_clip,
161 const int max_iter,
162 const cpl_boolean debug_data,
163 cpl_mask * bad_pixel_map,
164 cpl_frameset ** product_frames);
165
166static cpl_error_code
167eris_handle_cube_format_case(cpl_frameset * frameset,
168 const cpl_parameterlist * parlist,
169 cpl_frameset * raw_frames,
170 const char * recipe_name,
171 const char* combine_method,
172 const double sigma_clip,
173 const int max_iter,
174 const cpl_boolean debug_data,
175 cpl_mask * bad_pixel_map,
176 cpl_frameset ** product_frames);
177
178
179static cpl_error_code
180eris_handle_cube_format_robert(cpl_frameset * frameset,
181 const cpl_parameterlist * parlist,
182 cpl_frameset * raw_frames,
183 const char * recipe_name,
184 const char* combine_method,
185 const double sigma_clip,
186 const int max_iter,
187 const cpl_boolean debug_data,
188 cpl_mask * bad_pixel_map,
189 cpl_frameset ** product_frames);
190
191static located_imagelist *
192eris_create_final_skysub_products2(const hdrl_imagelist* science_iml, const cpl_frameset* raw_frames,
193 const hdrl_image* sky_flat, cpl_image* confidence, cpl_frameset** products);
194
195
196static located_imagelist *
197eris_create_final_skysub_products3(const hdrl_imagelist* science_iml, const cpl_frameset* frameset,
198 const cpl_frameset* raw_frames, const hdrl_image* sky_flat, cpl_image* confidence,
199 cpl_frameset** products);
200
201static located_imagelist *
202eris_create_final_skysub_products_robert(const cpl_frameset* frameset, const cpl_frameset* raw_frames, const hdrl_image* sky_flat,
203 cpl_image* confidence, cpl_frameset** products);
204
205static const int debug = 0;
206/*-----------------------------------------------------------------------------
207 Functions
208 -----------------------------------------------------------------------------*/
209
210/*----------------------------------------------------------------------------*/
226/*----------------------------------------------------------------------------*/
227cpl_error_code eris_nix_img_supersky_run(
228 cpl_frameset * frameset,
229 const cpl_parameterlist * parlist,
230 const char* recipe_name,
231 cpl_frameset ** product_frames)
232{
233 cpl_frameset * raw_frames = NULL;
234 const cpl_frame * bpm_frame = NULL;
235 cpl_mask * bad_pixel_map = NULL;
236
237 char * param_name = NULL;
238 const char* context = "eris.eris_nix_img_supersky";
239 const cpl_parameter* par;
240 int num_frames = 0;
241
242 cpl_error_code error = CPL_ERROR_NONE;
243
244 cpl_ensure_code(frameset, CPL_ERROR_NULL_INPUT);
245 cpl_ensure_code(parlist, CPL_ERROR_NULL_INPUT);
246 cpl_ensure_code(product_frames, CPL_ERROR_NULL_INPUT);
247
248 /* Extract parameters */
249 param_name = cpl_sprintf("%s.combine_method", context);
250 par = cpl_parameterlist_find_const(parlist, param_name);
251 const char* combine_method = cpl_parameter_get_string(par);
252 cpl_free(param_name);
253
254 param_name = cpl_sprintf("%s.sigma_clip", context);
255 par = cpl_parameterlist_find_const(parlist, param_name);
256 const double sigma_clip = cpl_parameter_get_double(par);
257 cpl_free(param_name);
258
259 param_name = cpl_sprintf("%s.max_iter", context);
260 par = cpl_parameterlist_find_const(parlist, param_name);
261 const int max_iter = cpl_parameter_get_int(par);
262 cpl_free(param_name);
263
264 param_name = cpl_sprintf("%s.save_skysub", context);
265 par = cpl_parameterlist_find_const(parlist, param_name);
266 const cpl_boolean save_skysub = cpl_parameter_get_bool(par);
267 cpl_free(param_name);
268
269 param_name = cpl_sprintf("%s.debug_data", context);
270 par = cpl_parameterlist_find_const(parlist, param_name);
271 const cpl_boolean debug_data = cpl_parameter_get_bool(par);
272 cpl_free(param_name);
273
274 /* Extract science frames */
275 raw_frames = eris_dfs_extract_frames_with_tag (frameset, ERIS_NIX_CAL_DET_OBJECT_JITTER_PRO_CATG);
276
277 if (!raw_frames || cpl_frameset_get_size(raw_frames) == 0) {
278 raw_frames = eris_dfs_extract_frames_with_tag (frameset, ERIS_NIX_CAL_DET_STD_JITTER_PRO_CATG);
279 }
280 if (!raw_frames || cpl_frameset_get_size(raw_frames) == 0) {
281 cpl_msg_error(cpl_func, "No science frames found with tag %s or %s",
282 ERIS_NIX_CAL_DET_OBJECT_JITTER_PRO_CATG, ERIS_NIX_CAL_DET_STD_JITTER_PRO_CATG);
283 error = CPL_ERROR_DATA_NOT_FOUND;
284 goto cleanup;
285 }
286 num_frames = cpl_frameset_get_size(raw_frames);
287 cpl_msg_info(cpl_func, "Processing %d science frames", num_frames);
288
289 /* Load bad pixel map */
290 bpm_frame = cpl_frameset_find_const(frameset, ERIS_NIX_IMG_SUPERSKY_BPM);
291 if (!bpm_frame) {
292 cpl_msg_warning(cpl_func, "No bad pixel map found with tag %s", ERIS_NIX_IMG_SUPERSKY_BPM);
293 }
294 cpl_image * bpm_image = NULL;
295 //cpl_mask *inp_mask = NULL;
296 if(bpm_frame) {
297 bpm_image = cpl_image_load(cpl_frame_get_filename(bpm_frame), CPL_TYPE_INT, 0, 0);
298
299 /* Convert BPM image to mask (0=good, >0=bad) */
300 bad_pixel_map = cpl_mask_threshold_image_create(bpm_image, 0.5, DBL_MAX);
301 if(debug) {
302 cpl_mask_save(bad_pixel_map, "bad_pixel_map.fits", NULL, CPL_IO_DEFAULT);
303 }
304 //inp_mask = cpl_mask_threshold_image_create(bpm_image, -0.5, 0.5);
305 }
306 cpl_boolean format_is_cube_of_same_size = CPL_TRUE;
307 cpl_boolean test_cubes_slices = CPL_TRUE;
308 format_is_cube_of_same_size = eris_check_format_is_cube_of_same_size(raw_frames);
309
310 //char* fname;
311 *product_frames = cpl_frameset_new();
312 if(!format_is_cube_of_same_size) {
313 num_frames = cpl_frameset_get_size(raw_frames);
314 cpl_msg_info(cpl_func, "Processing1 %d science frames", num_frames);
315
316 eris_handle_simple_format_case(frameset, parlist, raw_frames, recipe_name, combine_method,
317 sigma_clip, max_iter, debug_data, bad_pixel_map, product_frames);
318
319 } else {
320 /*
321 eris_handle_cube_format_case(frameset, parlist, raw_frames, recipe_name, combine_method,
322 sigma_clip, max_iter, debug_data, bad_pixel_map, product_frames);
323 */
324 eris_handle_cube_format_robert(frameset, parlist, raw_frames, recipe_name, combine_method,
325 sigma_clip, max_iter, debug_data, bad_pixel_map, product_frames);
326
327 }
328
329 eris_print_rec_status(700);
330
331 eris_print_rec_status(800);
332 cpl_msg_info(cpl_func, "Super-sky recipe completed successfully");
333
334cleanup:
335 eris_print_rec_status(600);
336 if (raw_frames) cpl_frameset_delete(raw_frames);
337 if (bad_pixel_map) cpl_mask_delete(bad_pixel_map);
338 //if (object_jitters) enu_located_imagelist_delete(object_jitters);
339 //if (sky_1st) enu_located_imagelist_delete(sky_1st);
340 //if(inp_mask) cpl_mask_delete(inp_mask);
341 //if(bpm_image !=NULL) cpl_image_delete(bpm_image);
342
343 eris_check_error_code(cpl_func);
344 return error;
345}
346
347
348
349static cpl_error_code
350eris_handle_simple_format_case(cpl_frameset * frameset,
351 const cpl_parameterlist * parlist,
352 cpl_frameset * raw_frames,
353 const char * recipe_name,
354 const char* combine_method,
355 const double sigma_clip,
356 const int max_iter,
357 const cpl_boolean debug_data,
358 cpl_mask * bad_pixel_map,
359 cpl_frameset ** product_frames){
360
361 hdrl_imagelist * science_images = NULL;
362 hdrl_image * sky_flat_0 = NULL;
363 hdrl_image * sky_flat_final = NULL;
364 cpl_mask ** source_masks = NULL;
365 hdrl_imagelist * skysub_images_0 = NULL;
366 /* Load science images as HDRL imagelist */
367 science_images = eris_crea_imagelist_all(raw_frames, bad_pixel_map);
368
369 /* Step 1: Create first-pass super-sky */
370 cpl_msg_info(cpl_func, "Creating first-pass super-sky flat...");
371 eris_print_rec_status(100);
372 sky_flat_0 = eris_create_sky_flat_first(science_images, bad_pixel_map, combine_method, sigma_clip,
373 max_iter, -1);
374
375 eris_print_rec_status(200);
376
377 /* Step 2: Subtract first-pass sky from each frame */
378 eris_first_sky_sub(science_images, sky_flat_0, &skysub_images_0);
379 if (sky_flat_0) hdrl_image_delete(sky_flat_0);
380
381 eris_print_rec_status(300);
382 /* Step 3: Detect sources in sky-subtracted frames */
383 cpl_size num_tot_images = hdrl_imagelist_get_size(skysub_images_0);
384 cpl_msg_info(cpl_func,"num_tot_images4: %lld",num_tot_images);
385 source_masks = cpl_calloc(num_tot_images, sizeof(cpl_mask *));
386 cpl_image * bpm_image = NULL;
387 char* prefix = cpl_sprintf("%s.catalogue",recipe_name);
388 eris_locate_and_mask_sources(skysub_images_0, bpm_image, parlist, prefix, source_masks);
389 cpl_free(prefix);
390 /* Step 4: Create final super-sky with source masks */
391 eris_print_rec_status(400);
392 cpl_msg_info(cpl_func, "Creating final super-sky flat with source masking...");
393 cpl_image* contrib = NULL;
394 sky_flat_final = eris_create_sky_flat_final(science_images, source_masks,
395 combine_method, sigma_clip, max_iter, -1, &contrib);
396
397 if (!sky_flat_final) {
398 cpl_msg_error(cpl_func, "Failed to create final sky flat");
399 return CPL_ERROR_ILLEGAL_OUTPUT;
400 }
401
402 if(sky_flat_final) {
403 /* Step 5: Save products */
404 eris_print_rec_status(500);
405
406 eris_save_sky_flat_final(raw_frames, sky_flat_final, product_frames);
407 }
408
409 /* Save sky-subtracted frames if requested */
410 //if (debug) {
411 //eris_create_final_skysub_products(science_images, sky_flat_final, product_frames);
412 located_imagelist * skysub = eris_create_final_skysub_products3(science_images, frameset, raw_frames, sky_flat_final, contrib, product_frames);
413 //char * preface = cpl_sprintf("%s_%s", "jitter", "skysub");
414 /* save the sky-subtracted jitter images to FITS files */
415 for (cpl_size i = 0; i < skysub->size; i++) {
416
417 // add necessary things to output header. RA and DEC are needed as
418 // would otherwise be overwritten by DFS
419 cpl_propertylist * applist = cpl_propertylist_new();
420 cpl_propertylist_update_string(applist, CPL_DFS_PRO_CATG, ERIS_NIX_SKYSUB_OBJECT_JITTER_PRO_CATG);
421 cpl_propertylist_update_string(applist, "PRODCATG",
422 "ANCILLARY.IMAGE");
423 cpl_propertylist_copy_property_regexp(applist,
424 skysub->limages[i]->plist,
425 "RA|DEC",
426 CPL_FALSE);
427
428 // set RA, DEC from WCS rather than let dfs copy it from ano†her frame
429 {
430 cpl_wcs * wcs = cpl_wcs_new_from_propertylist(
431 skysub->limages[i]->plist);
432 double ra = 0.0;
433 double dec = 0.0;
434 enu_get_ra_dec(wcs, &ra, &dec);
435
436 cpl_propertylist_update_double(applist, "RA", ra);
437 cpl_propertylist_update_double(applist, "DEC", dec);
438
439 cpl_wcs_delete(wcs);
440 }
441
442 // Generate output file name and write the file
443 char * out_fname = enu_repreface(cpl_frame_get_filename(
444 skysub->limages[i]->
445 frame),
446 "skysub");
447
448 enu_dfs_save_limage(frameset,
449 parlist,
450 raw_frames,
451 CPL_TRUE,
452 skysub->limages[i],
453 recipe_name,
454 skysub->limages[i]->frame,
455 applist,
456 PACKAGE "/" PACKAGE_VERSION,
457 out_fname);
458
459 cpl_free(out_fname);
460 cpl_propertylist_delete(applist);
461 }
462
463
464 /*
465 char * preface = cpl_sprintf("%s", "");
466 cpl_frameset_dump(frameset,stdout);
467 enu_debug_limlist_save(debug_data, skysub, preface, recipe_name, frameset, parlist,
468 raw_frames);
469 cpl_free(preface);
470 */
471 if (skysub) enu_located_imagelist_delete(skysub);
472 //}
473
474
475 if(contrib) cpl_image_delete(contrib);
476 if (science_images) hdrl_imagelist_delete(science_images);
477 if (sky_flat_final) hdrl_image_delete(sky_flat_final);
478 if (skysub_images_0) hdrl_imagelist_delete(skysub_images_0);
479 if (source_masks != NULL) {
480 for (int i = 0; i < num_tot_images; i++) {
481 if (source_masks[i]) cpl_mask_delete(source_masks[i]);
482 }
483 cpl_free(source_masks);
484 }
485
486 eris_check_error_code(cpl_func);
487 return cpl_error_get_code();
488}
489
490
491
492static cpl_error_code
493eris_handle_simple_format_case2(cpl_frameset * frameset,
494 const cpl_parameterlist * parlist,
495 const char * recipe_name,
496 const char* combine_method,
497 const double sigma_clip,
498 const int max_iter,
499 const cpl_boolean debug_data,
500 cpl_mask * bad_pixel_map,
501 cpl_frameset ** product_frames){
502
503 hdrl_imagelist * science_images = NULL;
504 hdrl_image * sky_flat_0 = NULL;
505 hdrl_image * sky_flat_final = NULL;
506 located_imagelist * object_jitters = NULL;
507 located_imagelist * sky_1st = NULL;
508 cpl_mask ** source_masks = NULL;
509 hdrl_imagelist * skysub_images_0 = NULL;
510 cpl_frameset * raw_frames = cpl_frameset_new();
511 /* Load science images as HDRL imagelist */
512
513
514 object_jitters = enu_limlist_load_from_frameset(frameset,
515 ERIS_NIX_CAL_DET_OBJECT_JITTER_PRO_CATG, raw_frames);
516 cpl_msg_info(cpl_func, "Processing3 %lld science frames", cpl_frameset_get_size(raw_frames));
517 /* Extract science frames */
518
519 if (!object_jitters) {
520 object_jitters = enu_limlist_load_from_frameset(frameset,
521 ERIS_NIX_CAL_DET_STD_JITTER_PRO_CATG, raw_frames);
522 }
523 if (!raw_frames || cpl_frameset_get_size(raw_frames) == 0) {
524 cpl_msg_error(cpl_func, "No science frames found with tag %s or %s",
525 ERIS_NIX_CAL_DET_OBJECT_JITTER_PRO_CATG, ERIS_NIX_CAL_DET_STD_JITTER_PRO_CATG);
526
527 return CPL_ERROR_DATA_NOT_FOUND;
528 }
529 cpl_msg_info(cpl_func,"n0: %lld",cpl_frameset_get_size(raw_frames));
530 science_images = eris_crea_imagelist_all(raw_frames, bad_pixel_map);
531 cpl_msg_info(cpl_func,"n1: %lld",hdrl_imagelist_get_size(science_images));
532 /* Step 1: Create first-pass super-sky */
533 cpl_msg_info(cpl_func, "Creating first-pass super-sky flat...");
534 /*
535 sky_flat_0 = eris_create_sky_flat_first(science_images, bad_pixel_map, combine_method, sigma_clip,
536 max_iter, -1);
537 */
538 eris_print_rec_status(100);
539 /*
540 if(debug) {
541 char * preface = cpl_sprintf("%s_%s", "jitters", "input");
542 enu_debug_limlist_save(debug_data, object_jitters, preface, recipe_name, frameset, parlist,
543 raw_frames);
544 cpl_free(preface);
545 }
546 */
547
548
549 sky_flat_0 = eris_create_sky_flat_first2(object_jitters, &sky_1st, bad_pixel_map, combine_method, sigma_clip,
550 max_iter, -1);
551 eris_print_rec_status(150);
552 cpl_msg_info(cpl_func,"n2: %lld",object_jitters->size);
553
554
555 if(debug) {
556 char * preface = cpl_sprintf("%s_%s", "jitters", "sky_1st");
557 enu_debug_limlist_save(debug_data, sky_1st, preface, recipe_name, frameset, parlist,
558 raw_frames);
559 cpl_free(preface);
560 }
561
562
563 eris_print_rec_status(170);
564
565 /* Step 2: Subtract first-pass sky from each frame */
566 eris_first_sky_sub(science_images, sky_flat_0, &skysub_images_0);
567 if (sky_flat_0) hdrl_image_delete(sky_flat_0);
568 cpl_msg_info(cpl_func,"n3: %lld",hdrl_imagelist_get_size(science_images));
569 cpl_msg_info(cpl_func,"n4: %lld",hdrl_imagelist_get_size(skysub_images_0));
570
571 eris_print_rec_status(300);
572 /* Step 3: Detect sources in sky-subtracted frames */
573 cpl_size num_tot_images = hdrl_imagelist_get_size(skysub_images_0);
574 cpl_msg_info(cpl_func,"num_tot_images1: %lld",num_tot_images);
575 source_masks = cpl_calloc(num_tot_images, sizeof(cpl_mask *));
576 cpl_image * bpm_image = NULL;
577 char* prefix = cpl_sprintf("%s.catalogue",recipe_name);
578 //eris_locate_and_mask_sources(skysub_images_0, bpm_image, parlist, recipe_name, source_masks);
579 eris_locate_and_mask_sources2(sky_1st, bpm_image, parlist, prefix, source_masks);
580 cpl_free(prefix);
581 eris_print_rec_status(400);
582 /* Step 4: Create final super-sky with source masks */
583 cpl_msg_info(cpl_func, "Creating final super-sky flat with source masking...");
584 //sky_flat_final = eris_create_sky_flat_final(science_images, source_masks,
585 // combine_method, sigma_clip, max_iter, -1);
586 if (sky_1st) enu_located_imagelist_delete(sky_1st);
587 located_imagelist* sky_2nd = enu_located_imagelist_duplicate(object_jitters);
588 sky_flat_final = eris_create_sky_flat_final2(object_jitters, &sky_2nd, source_masks, combine_method, sigma_clip, max_iter, -1);
589
590 eris_print_rec_status(450);
591
592 eris_print_rec_status(460);
593
594 if (!sky_flat_final) {
595 cpl_msg_error(cpl_func, "Failed to create final sky flat");
596 return CPL_ERROR_ILLEGAL_OUTPUT;
597 }
598 /*
599 if(debug) {
600 char * preface = cpl_sprintf("%s_%s", "jitters", "sky_2nd");
601 enu_debug_limlist_save(debug_data, sky_2nd, preface, recipe_name, frameset, parlist,
602 raw_frames);
603 cpl_free(preface);
604 }
605 */
606
607 if(sky_flat_final) {
608 /* Step 5: Save products */
609 eris_print_rec_status(500);
610
611 eris_save_sky_flat_final(raw_frames, sky_flat_final, product_frames);
612 }
613
614 /* Save sky-subtracted frames if requested */
615 if (debug) {
616 eris_create_final_skysub_products(science_images, sky_flat_final, product_frames);
617 }
618
619 if (science_images) hdrl_imagelist_delete(science_images);
620 if (sky_flat_final) hdrl_image_delete(sky_flat_final);
621 if (skysub_images_0) hdrl_imagelist_delete(skysub_images_0);
622 //if (object_jitters) enu_located_imagelist_delete(object_jitters);
623 //if (sky_2nd) enu_located_imagelist_delete(sky_2nd);
624 if (source_masks != NULL) {
625 for (int i = 0; i < num_tot_images; i++) {
626 if (source_masks[i]) cpl_mask_delete(source_masks[i]);
627 }
628 cpl_free(source_masks);
629 }
630 cpl_frameset_delete(raw_frames);
631 eris_check_error_code(cpl_func);
632 return cpl_error_get_code();
633}
634
635
636static cpl_error_code
637eris_handle_cube_format_case(cpl_frameset * frameset,
638 const cpl_parameterlist * parlist,
639 cpl_frameset * raw_frames,
640 const char * recipe_name,
641 const char* combine_method,
642 const double sigma_clip,
643 const int max_iter,
644 const cpl_boolean debug_data,
645 cpl_mask * bad_pixel_map,
646 cpl_frameset ** product_frames)
647{
648
649 hdrl_imagelist * science_images = NULL;
650 hdrl_image * sky_flat_0 = NULL;
651 hdrl_image * sky_flat_final = NULL;
652 cpl_mask ** source_masks = NULL;
653 hdrl_imagelist * skysub_images_0 = NULL;
654 cpl_msg_warning(cpl_func,"cube format, process slice by slice");
655
656 const cpl_frame * frame = cpl_frameset_get_position_const(raw_frames, 0);
657 const char * filename = cpl_frame_get_filename(frame);
658 cpl_imagelist * data_iml = cpl_imagelist_load(filename, CPL_TYPE_DOUBLE, 1);
659 cpl_size nplanes = cpl_imagelist_get_size(data_iml);
660 cpl_size num_frames = cpl_frameset_get_size(raw_frames);
661 for(cpl_size k = 0; k < nplanes; k++) {
662 science_images = hdrl_imagelist_new();
663 for(cpl_size kkk = 0; kkk < num_frames; kkk++) {
664 frame = cpl_frameset_get_position_const(raw_frames, kkk);
665 filename = cpl_frame_get_filename(frame);
666 eris_load_and_error_cube_slice(filename, bad_pixel_map, k, kkk, &science_images);
667 }
668
669 /* Step 1: Create first-pass super-sky */
670 cpl_msg_info(cpl_func, "Creating first-pass super-sky flat...");
671 sky_flat_0 = eris_create_sky_flat_first(science_images, bad_pixel_map, combine_method, sigma_clip, max_iter, k);
672
673 /* Step 2: Subtract first-pass sky from each frame */
674 eris_first_sky_sub(science_images, sky_flat_0, &skysub_images_0);
675
676 eris_print_rec_status(300);
677 /* Step 3: Detect sources in sky-subtracted frames */
678 cpl_size num_tot_images = hdrl_imagelist_get_size(skysub_images_0);
679 cpl_msg_info(cpl_func,"num_tot_images2: %lld",num_tot_images);
680 source_masks = cpl_calloc(num_tot_images, sizeof(cpl_mask *));
681 cpl_image * bpm_image = NULL;
682 char* prefix = cpl_sprintf("%s.catalogue",recipe_name);
683 eris_locate_and_mask_sources(skysub_images_0, bpm_image, parlist, prefix, source_masks);
684 cpl_free(prefix);
685 eris_print_rec_status(400);
686 /* Step 4: Create final super-sky with source masks */
687 cpl_msg_info(cpl_func, "Creating final super-sky flat with source masking...");
688 cpl_image* contrib = NULL;
689 sky_flat_final = eris_create_sky_flat_final(science_images, source_masks,
690 combine_method, sigma_clip, max_iter, k, &contrib);
691
692 //if (debug) {
693 eris_create_final_skysub_products(science_images, sky_flat_final, product_frames);
694 //}
695
696 if (sky_flat_0) hdrl_image_delete(sky_flat_0);
697 if (science_images) hdrl_imagelist_delete(science_images);
698 if (sky_flat_final) hdrl_image_delete(sky_flat_final);
699 if (skysub_images_0) hdrl_imagelist_delete(skysub_images_0);
700 if (contrib) cpl_image_delete(contrib);
701 if (source_masks != NULL) {
702 for (int i = 0; i < num_tot_images; i++) {
703 if (source_masks[i]) cpl_mask_delete(source_masks[i]);
704 }
705 cpl_free(source_masks);
706 }
707 }
708 cpl_imagelist_delete(data_iml);
709
710 eris_check_error_code(cpl_func);
711 return cpl_error_get_code();
712
713}
714
715static cpl_error_code
716eris_handle_cube_format_robert(cpl_frameset * frameset,
717 const cpl_parameterlist * parlist,
718 cpl_frameset * raw_frames,
719 const char * recipe_name,
720 const char* combine_method,
721 const double sigma_clip,
722 const int max_iter,
723 const cpl_boolean debug_data,
724 cpl_mask * bad_pixel_map,
725 cpl_frameset ** product_frames)
726{
727
728 hdrl_imagelist * science_images = NULL;
729 hdrl_image * sky_flat_0 = NULL;
730 hdrl_image * sky_flat_final = NULL;
731 cpl_mask ** source_masks = NULL;
732 hdrl_imagelist * skysub_images_0 = NULL;
733 cpl_msg_warning(cpl_func,"cube format, collapse cubes and combine images");
734
735 const cpl_frame * frame = cpl_frameset_get_position_const(raw_frames, 0);
736 const char * filename = cpl_frame_get_filename(frame);
737
738
739 cpl_size num_frames = cpl_frameset_get_size(raw_frames);
740
741
742 science_images = hdrl_imagelist_new();
743 for(cpl_size kkk = 0; kkk < num_frames; kkk++) {
744 frame = cpl_frameset_get_position_const(raw_frames, kkk);
745 filename = cpl_frame_get_filename(frame);
746
747 eris_load_and_error_cube_robert(filename, bad_pixel_map, combine_method, sigma_clip, max_iter, kkk,
748 &science_images);
749
750 }
751
752 /* Step 1: Create first-pass super-sky */
753 cpl_msg_info(cpl_func, "Creating first-pass super-sky flat...");
754 sky_flat_0 = eris_create_sky_flat_first(science_images, bad_pixel_map, combine_method, sigma_clip, max_iter, 0);
755
756 /* Step 2: Subtract first-pass sky from each frame */
757 eris_first_sky_sub(science_images, sky_flat_0, &skysub_images_0);
758 if (sky_flat_0) hdrl_image_delete(sky_flat_0);
759
760 eris_print_rec_status(300);
761 /* Step 3: Detect sources in sky-subtracted frames */
762 cpl_size num_tot_images = hdrl_imagelist_get_size(skysub_images_0);
763 cpl_msg_info(cpl_func,"num_tot_images2: %lld",num_tot_images);
764 source_masks = cpl_calloc(num_tot_images, sizeof(cpl_mask *));
765 cpl_image * bpm_image = NULL;
766 char* prefix = cpl_sprintf("%s.catalogue",recipe_name);
767 eris_locate_and_mask_sources(skysub_images_0, bpm_image, parlist, prefix, source_masks);
768 hdrl_imagelist_delete(skysub_images_0);
769 cpl_free(prefix);
770
771 eris_print_rec_status(400);
772 /* Step 4: Create final super-sky with source masks */
773 cpl_msg_info(cpl_func, "Creating final super-sky flat with source masking...");
774 cpl_image* contrib = NULL;
775 sky_flat_final = eris_create_sky_flat_final(science_images, source_masks,
776 combine_method, sigma_clip, max_iter, 0, &contrib);
777
778 if (science_images) hdrl_imagelist_delete(science_images);
779 if (source_masks != NULL) {
780 for (int i = 0; i < num_tot_images; i++) {
781 if (source_masks[i]) cpl_mask_delete(source_masks[i]);
782 }
783 cpl_free(source_masks);
784 }
785
786 eris_create_final_skysub_products_robert(frameset, raw_frames, sky_flat_final, contrib, product_frames);
787 //eris_create_final_skysub_products(science_images, sky_flat_final, product_frames);
788 //cpl_msg_info(cpl_func,"about to exit");
789
790 if (sky_flat_final) hdrl_image_delete(sky_flat_final);
791 if (contrib) cpl_image_delete(contrib);
792
793 eris_check_error_code(cpl_func);
794 return cpl_error_get_code();
795
796}
797
798static cpl_error_code
799eris_handle_cube_format_case2(cpl_frameset * frameset,
800 const cpl_parameterlist * parlist,
801 cpl_frameset * raw_frames,
802 const char * recipe_name,
803 const char* combine_method,
804 const double sigma_clip,
805 const int max_iter,
806 const cpl_boolean debug_data,
807 cpl_mask * bad_pixel_map,
808 cpl_frameset ** product_frames)
809{
810
811 hdrl_imagelist * science_images = NULL;
812 hdrl_image * sky_flat_0 = NULL;
813 hdrl_image * sky_flat_final = NULL;
814 cpl_mask ** source_masks = NULL;
815 hdrl_imagelist * skysub_images_0 = NULL;
816 cpl_msg_warning(cpl_func,"cube format, process slice by slice");
817
818 const cpl_frame * frame = cpl_frameset_get_position_const(raw_frames, 0);
819 const char * filename = cpl_frame_get_filename(frame);
820 cpl_imagelist * data_iml = cpl_imagelist_load(filename, CPL_TYPE_DOUBLE, 1);
821 cpl_size nplanes = cpl_imagelist_get_size(data_iml);
822 cpl_size num_frames = cpl_frameset_get_size(raw_frames);
823 for(cpl_size k = 0; k < nplanes; k++) {
824 science_images = hdrl_imagelist_new();
825 for(cpl_size kkk = 0; kkk < num_frames; kkk++) {
826 frame = cpl_frameset_get_position_const(raw_frames, kkk);
827 filename = cpl_frame_get_filename(frame);
828 eris_load_and_error_cube_slice(filename, bad_pixel_map, k, kkk, &science_images);
829 }
830
831 /* Step 1: Create first-pass super-sky */
832 cpl_msg_info(cpl_func, "Creating first-pass super-sky flat...");
833 sky_flat_0 = eris_create_sky_flat_first(science_images, bad_pixel_map, combine_method, sigma_clip, max_iter, k);
834
835 /* Step 2: Subtract first-pass sky from each frame */
836 eris_first_sky_sub(science_images, sky_flat_0, &skysub_images_0);
837
838 eris_print_rec_status(300);
839 /* Step 3: Detect sources in sky-subtracted frames */
840 cpl_size num_tot_images = hdrl_imagelist_get_size(skysub_images_0);
841 cpl_msg_info(cpl_func,"num_tot_images2: %lld",num_tot_images);
842 source_masks = cpl_calloc(num_tot_images, sizeof(cpl_mask *));
843 cpl_image * bpm_image = NULL;
844 char* prefix = cpl_sprintf("%s.catalogue",recipe_name);
845 eris_locate_and_mask_sources(skysub_images_0, bpm_image, parlist, prefix, source_masks);
846 cpl_free(prefix);
847 eris_print_rec_status(400);
848 /* Step 4: Create final super-sky with source masks */
849 cpl_msg_info(cpl_func, "Creating final super-sky flat with source masking...");
850 cpl_image* contrib = NULL;
851 sky_flat_final = eris_create_sky_flat_final(science_images, source_masks,
852 combine_method, sigma_clip, max_iter, k, &contrib);
853
854 //if (debug) {
855 eris_create_final_skysub_products(science_images, sky_flat_final, product_frames);
856 //}
857
858 if (sky_flat_0) hdrl_image_delete(sky_flat_0);
859 if (science_images) hdrl_imagelist_delete(science_images);
860 if (sky_flat_final) hdrl_image_delete(sky_flat_final);
861 if (skysub_images_0) hdrl_imagelist_delete(skysub_images_0);
862 if (contrib) cpl_image_delete(contrib);
863 if (source_masks != NULL) {
864 for (int i = 0; i < num_tot_images; i++) {
865 if (source_masks[i]) cpl_mask_delete(source_masks[i]);
866 }
867 cpl_free(source_masks);
868 }
869 }
870 cpl_imagelist_delete(data_iml);
871
872 eris_check_error_code(cpl_func);
873 return cpl_error_get_code();
874
875}
876
877
878
879static cpl_error_code
880eris_save_sky_flat_final( const cpl_frameset* raw_frames, const hdrl_image* sky_flat_final, cpl_frameset** product_frames)
881{
882
883 /* Step 5: Save products */
884 cpl_size num_frames = cpl_frameset_get_size(raw_frames);
885 cpl_msg_info(cpl_func,"num_frames:%lld",num_frames);
886 eris_print_rec_status(500);
887 //*product_frames = cpl_frameset_new();
888 hdrl_value med_sky_final = hdrl_image_get_median(sky_flat_final);
889
890 /* Save master sky flat */
891 cpl_propertylist* plist = cpl_propertylist_load(cpl_frame_get_filename(
892 cpl_frameset_get_position_const(raw_frames, 0)), 0);
893 cpl_propertylist_update_string(plist, CPL_DFS_PRO_CATG, ERIS_NIX_IMG_SUPERSKY_SKYFLAT);
894 cpl_propertylist_update_double(plist, "ESO QC SKY MEDIAN", med_sky_final.data);
895 cpl_propertylist_update_int(plist, "ESO QC NFRAMES", num_frames);
896 eris_print_rec_status(600);
897 char* fname = cpl_sprintf("sky_flat_2nd.fits");
898
899 cpl_image_save(hdrl_image_get_image_const(sky_flat_final), fname,CPL_TYPE_DOUBLE, plist, CPL_IO_CREATE);
900 cpl_image_save(hdrl_image_get_error_const(sky_flat_final), fname,CPL_TYPE_DOUBLE, plist, CPL_IO_EXTEND);
901 cpl_mask_save(hdrl_image_get_mask_const(sky_flat_final), fname, plist, CPL_IO_EXTEND);
902
903 cpl_frame * product = cpl_frame_new();
904 cpl_frame_set_filename(product, fname);
905 cpl_frame_set_tag(product, ERIS_NIX_IMG_SUPERSKY_SKYFLAT);
906 cpl_frame_set_type(product, CPL_FRAME_TYPE_IMAGE);
907 cpl_frame_set_group(product, CPL_FRAME_GROUP_PRODUCT);
908 cpl_frameset_insert(*product_frames, product);
909
910 cpl_free(fname);
911 cpl_propertylist_delete(plist);
912
913 eris_check_error_code(cpl_func);
914 return cpl_error_get_code();
915
916}
917
918static cpl_error_code
919eris_save_sky_flat_final2( const cpl_frameset* raw_frames, const hdrl_image* sky_flat_final, cpl_frameset** product_frames)
920{
921
922 /* Step 5: Save products */
923 cpl_size num_frames = cpl_frameset_get_size(raw_frames);
924 eris_print_rec_status(500);
925 //*product_frames = cpl_frameset_new();
926 hdrl_value med_sky_final = hdrl_image_get_median(sky_flat_final);
927
928 /* Save master sky flat */
929 cpl_propertylist* plist = cpl_propertylist_load(cpl_frame_get_filename(
930 cpl_frameset_get_position_const(raw_frames, 0)), 0);
931 cpl_propertylist_update_string(plist, CPL_DFS_PRO_CATG, ERIS_NIX_IMG_SUPERSKY_SKYFLAT);
932 cpl_propertylist_update_double(plist, "ESO QC SKY MEDIAN", med_sky_final.data);
933 cpl_propertylist_update_int(plist, "ESO QC NFRAMES", num_frames);
934 eris_print_rec_status(600);
935 char* fname = cpl_sprintf("sky_flat_2nd.fits");
936
937 cpl_image_save(hdrl_image_get_image_const(sky_flat_final), fname,CPL_TYPE_DOUBLE, plist, CPL_IO_CREATE);
938 cpl_image_save(hdrl_image_get_error_const(sky_flat_final), fname,CPL_TYPE_DOUBLE, plist, CPL_IO_EXTEND);
939 cpl_mask_save(hdrl_image_get_mask_const(sky_flat_final), fname, plist, CPL_IO_EXTEND);
940
941 cpl_frame * product = cpl_frame_new();
942 cpl_frame_set_filename(product, fname);
943 cpl_frame_set_tag(product, ERIS_NIX_IMG_SUPERSKY_SKYFLAT);
944 cpl_frame_set_type(product, CPL_FRAME_TYPE_IMAGE);
945 cpl_frame_set_group(product, CPL_FRAME_GROUP_PRODUCT);
946 cpl_frameset_insert(*product_frames, product);
947
948 cpl_free(fname);
949 cpl_propertylist_delete(plist);
950
951 eris_check_error_code(cpl_func);
952 return cpl_error_get_code();
953
954}
955
956static hdrl_imagelist*
957eris_crea_imagelist_all(cpl_frameset* raw_frames, cpl_mask* bad_pixel_map)
958{
959
960 /* Load science images as HDRL imagelist */
961 hdrl_imagelist* science_images = hdrl_imagelist_new();
962 cpl_size num_frames = cpl_frameset_get_size(raw_frames);
963 cpl_size next = 0;
964 const char* format = NULL;
965 cpl_propertylist* pheader = NULL;
966 eris_print_rec_status(100);
967 cpl_size kkk = 0;
968 for (cpl_size kk = 0; kk < num_frames; kk++) {
969 const cpl_frame * frame = cpl_frameset_get_position_const(raw_frames, kk);
970 const char * filename = cpl_frame_get_filename(frame);
971 next = cpl_frame_get_nextensions(frame);
972 pheader = cpl_propertylist_load(filename,0);
973 format = eris_pfits_get_frame_format(pheader);
974
975 if (next < 4) {
976 eris_load_data_and_crea_error(filename, kk, bad_pixel_map, &science_images);
977 } else {
978 cpl_msg_info(cpl_func,"sci data with extentions %s", filename);
979 if(!strcmp(format,"single")) {
980 eris_load_data_and_error_simple(filename, kk, bad_pixel_map, &science_images);
981 } else {
982 eris_load_and_error_cube(filename, bad_pixel_map, &kkk, &science_images);
983 }/* end cube format case */
984 } /* end case with extensions */
985 cpl_propertylist_delete(pheader);
986 } /* end loop over input frames */
987
988 cpl_size nx = hdrl_image_get_size_x(hdrl_imagelist_get_const(science_images, 0));
989 cpl_size ny = hdrl_image_get_size_y(hdrl_imagelist_get_const(science_images, 0));
990 cpl_msg_info(cpl_func, "Image dimensions: %lld x %lld", nx, ny);
991
992
993 eris_check_error_code(cpl_func);
994 return science_images;
995}
996
997static cpl_boolean
998eris_check_format_is_cube_of_same_size(const cpl_frameset * raw_frames)
999{
1000
1001 cpl_boolean result = CPL_TRUE;
1002 cpl_size nplanes = 0;
1003 cpl_size nplanes_tmp = 0;
1004
1005 cpl_size num_frames = cpl_frameset_get_size(raw_frames);
1006 cpl_size next = 0;
1007 for (int kk = 0; kk < num_frames; kk++) {
1008 const cpl_frame * frame = cpl_frameset_get_position_const(raw_frames, kk);
1009 next = cpl_frame_get_nextensions(frame);
1010
1011 if(next > 3) {
1012 const char * fname = cpl_frame_get_filename(frame);
1013
1014 cpl_propertylist* xhead = cpl_propertylist_load(fname,1);
1015 if(cpl_propertylist_has(xhead,"CD3_3")) {
1016 cpl_imagelist * data_iml = cpl_imagelist_load(fname, CPL_TYPE_DOUBLE, 1);
1017
1018 nplanes_tmp = cpl_imagelist_get_size(data_iml);
1019 if(kk == 0) {
1020 nplanes = nplanes_tmp;
1021 } else if (nplanes_tmp != nplanes) {
1022 cpl_msg_info(cpl_func,"cubes not of same size");
1023 result = CPL_FALSE;
1024 }
1025 cpl_imagelist_delete(data_iml);
1026 } else {
1027 cpl_msg_info(cpl_func,"Input Data are simple format (images)");
1028 cpl_propertylist_delete(xhead);
1029 return CPL_FALSE;
1030 }
1031 cpl_propertylist_delete(xhead);
1032 } else {
1033 cpl_msg_info(cpl_func,"Data with less than 4 extensions");
1034 return CPL_FALSE;
1035 }
1036 }
1037 eris_check_error_code(cpl_func);
1038 return result;
1039}
1040
1041/*----------------------------------------------------------------------------*/
1050/*----------------------------------------------------------------------------*/
1051
1052static cpl_error_code
1053eris_create_final_skysub_products(const hdrl_imagelist* science_iml, const hdrl_image* sky_flat, cpl_frameset** products)
1054{
1055
1056 cpl_msg_info(cpl_func, "Saving sky-subtracted frames...");
1057 cpl_size num_tot_images = hdrl_imagelist_get_size(science_iml);
1058 hdrl_value med_sky = hdrl_image_get_median(sky_flat);
1059
1060 for (int i = 0; i < num_tot_images; i++) {
1061
1062 const hdrl_image * himg = hdrl_imagelist_get_const(science_iml, i);
1063 hdrl_value med_img = hdrl_image_get_median(himg);
1064 hdrl_value scale = {0, 0};
1065 scale.data = med_img.data / med_sky.data;
1066
1067 hdrl_image * skysub = hdrl_image_duplicate(himg);
1068 hdrl_image * scaled_sky = hdrl_image_duplicate(sky_flat);
1069 hdrl_image_mul_scalar(scaled_sky, scale);
1070 hdrl_image_sub_image(skysub, scaled_sky);
1071 hdrl_value med_skyub = hdrl_image_get_median(skysub);
1072 /* Create output filename */
1073 char outfile[256];
1074 snprintf(outfile, 256, "sci_skysub_2nd_%03d.fits", i);
1075
1076 /* Load and update header */
1077 cpl_msg_info(cpl_func,"ESO QC SKY_MED: %g",med_img.data);
1078 cpl_msg_info(cpl_func,"ESO QC SKYSUB_MED: %g", med_skyub.data);
1079 /*
1080 const cpl_frame * frame = cpl_frameset_get_position_const(raw_frames, i);
1081 plist = cpl_propertylist_load(cpl_frame_get_filename(frame), 0);
1082 cpl_propertylist_update_string(plist, CPL_DFS_PRO_CATG,
1083 ERIS_NIX_IMG_SUPERSKY_SKYSUB);
1084
1085 //cpl_propertylist_update_double(plist, "ESO QC SKY_MED", med_img.data);
1086 //cpl_propertylist_update_double(plist, "ESO QC SKYSUB_MED", med_skyub.data);
1087
1088 */
1089
1090 cpl_image_save(hdrl_image_get_image(skysub), outfile, CPL_TYPE_DOUBLE, NULL, CPL_IO_CREATE);
1091 cpl_image_save(hdrl_image_get_error(skysub), outfile, CPL_TYPE_DOUBLE, NULL, CPL_IO_EXTEND);
1092 cpl_mask_save(hdrl_image_get_mask(skysub), outfile, NULL, CPL_IO_EXTEND);
1093
1094 cpl_frame * product = cpl_frame_new();
1095 cpl_frame_set_filename(product, outfile);
1096 cpl_frame_set_tag(product, ERIS_NIX_IMG_SUPERSKY_SKYSUB);
1097 cpl_frame_set_type(product, CPL_FRAME_TYPE_IMAGE);
1098 cpl_frame_set_group(product, CPL_FRAME_GROUP_PRODUCT);
1099 cpl_frameset_insert(*products, product);
1100 //cpl_propertylist_delete(plist);
1101 hdrl_image_delete(skysub);
1102 hdrl_image_delete(scaled_sky);
1103 }
1104
1105 eris_check_error_code(cpl_func);
1106 return cpl_error_get_code();
1107}
1108
1109
1110
1111
1112
1113static located_imagelist *
1114eris_create_final_skysub_products2(const hdrl_imagelist* science_iml, const cpl_frameset* raw_frames, const hdrl_image* sky_flat,
1115 cpl_image* confidence, cpl_frameset** products)
1116{
1117
1118 cpl_msg_info(cpl_func, "Saving sky-subtracted frames...");
1119 cpl_size num_tot_images = hdrl_imagelist_get_size(science_iml);
1120 hdrl_value med_sky = hdrl_image_get_median(sky_flat);
1121 located_imagelist * result = enu_located_imagelist_new(num_tot_images);
1122 for (int i = 0; i < num_tot_images; i++) {
1123
1124 const hdrl_image * himg = hdrl_imagelist_get_const(science_iml, i);
1125 hdrl_value med_img = hdrl_image_get_median(himg);
1126 hdrl_value scale = {0, 0};
1127 scale.data = med_img.data / med_sky.data;
1128
1129 hdrl_image * skysub = hdrl_image_duplicate(himg);
1130 hdrl_image * scaled_sky = hdrl_image_duplicate(sky_flat);
1131 hdrl_image_mul_scalar(scaled_sky, scale);
1132 hdrl_image_sub_image(skysub, scaled_sky);
1133 hdrl_value med_skyub = hdrl_image_get_median(skysub);
1134 /* Create output filename */
1135 char outfile[256];
1136 snprintf(outfile, 256, "ERIS.2024.sci_skysub_2nd_%03d.fits", i);
1137
1138 /* Load and update header */
1139 cpl_msg_info(cpl_func,"ESO QC SKY_MED: %g",med_img.data);
1140 cpl_msg_info(cpl_func,"ESO QC SKYSUB_MED: %g", med_skyub.data);
1141 /*
1142 const cpl_frame * frame = cpl_frameset_get_position_const(raw_frames, i);
1143 plist = cpl_propertylist_load(cpl_frame_get_filename(frame), 0);
1144 cpl_propertylist_update_string(plist, CPL_DFS_PRO_CATG,
1145 ERIS_NIX_IMG_SUPERSKY_SKYSUB);
1146
1147 //cpl_propertylist_update_double(plist, "ESO QC SKY_MED", med_img.data);
1148 //cpl_propertylist_update_double(plist, "ESO QC SKYSUB_MED", med_skyub.data);
1149
1150 */
1151 cpl_frame * frame = cpl_frameset_get_position(raw_frames, i);
1152 cpl_propertylist* plist = cpl_propertylist_load(cpl_frame_get_filename(frame), 0);
1153 cpl_propertylist* dhead = cpl_propertylist_load(cpl_frame_get_filename(frame), 1);
1154 cpl_propertylist* ehead = cpl_propertylist_load(cpl_frame_get_filename(frame), 2);
1155 cpl_propertylist* qhead = cpl_propertylist_load(cpl_frame_get_filename(frame), 3);
1156 cpl_propertylist_save(plist, outfile, CPL_IO_CREATE);
1157 cpl_image_save(hdrl_image_get_image(skysub), outfile, CPL_TYPE_DOUBLE, dhead, CPL_IO_EXTEND);
1158 cpl_image_save(hdrl_image_get_error(skysub), outfile, CPL_TYPE_DOUBLE, ehead, CPL_IO_EXTEND);
1159 cpl_mask_save(hdrl_image_get_mask(skysub), outfile, qhead, CPL_IO_EXTEND);
1160 cpl_propertylist_delete(dhead);
1161 cpl_propertylist_delete(ehead);
1162 cpl_propertylist_delete(qhead);
1163
1164
1165
1166
1167 cpl_frame * product = cpl_frame_new();
1168 cpl_frame_set_filename(product, outfile);
1169 cpl_frame_set_tag(product, ERIS_NIX_IMG_SUPERSKY_SKYSUB);
1170 cpl_frame_set_type(product, CPL_FRAME_TYPE_IMAGE);
1171 cpl_frame_set_group(product, CPL_FRAME_GROUP_CALIB);
1172 cpl_frameset_insert(*products, product);
1173 located_image * limage = enu_load_limage_from_frame(product, &confidence, CPL_FALSE);
1174 limage->bkg = hdrl_image_duplicate(scaled_sky);
1175 limage->frame = product;
1176 //limage->confidence = cpl_image_duplicate(confidence);
1177 //limage->plist = cpl_propertylist_duplicate(plist);
1178 enu_located_imagelist_insert(result,limage, i);
1179
1180 cpl_propertylist_delete(plist);
1181 hdrl_image_delete(skysub);
1182 hdrl_image_delete(scaled_sky);
1183 }
1184
1185 eris_check_error_code(cpl_func);
1186 return result;
1187}
1188
1189static located_imagelist *
1190eris_create_final_skysub_products3(const hdrl_imagelist* science_iml, const cpl_frameset* frameset,
1191 const cpl_frameset* raw_frames, const hdrl_image* sky_flat,
1192 cpl_image* confidence, cpl_frameset** products)
1193{
1194
1195 cpl_msg_info(cpl_func, "Saving sky-subtracted frames...");
1196 cpl_size num_tot_images = hdrl_imagelist_get_size(science_iml);
1197 hdrl_value med_sky = hdrl_image_get_median(sky_flat);
1198 located_imagelist* object_jitters = NULL;
1199 //cpl_frameset_dump(frameset,stdout);
1200 cpl_size sx = hdrl_image_get_size_x(sky_flat);
1201 cpl_size sy = hdrl_image_get_size_y(sky_flat);
1202 if(cpl_frameset_count_tags(frameset, ERIS_NIX_CAL_DET_OBJECT_JITTER_PRO_CATG) > 0) {
1203 object_jitters = enu_limlist_load_from_frameset(frameset,
1204 ERIS_NIX_CAL_DET_OBJECT_JITTER_PRO_CATG, raw_frames);
1205 } else if(cpl_frameset_count_tags(frameset, ERIS_NIX_CAL_DET_STD_JITTER_PRO_CATG) > 0) {
1206 object_jitters = enu_limlist_load_from_frameset(frameset,
1207 ERIS_NIX_CAL_DET_STD_JITTER_PRO_CATG, raw_frames);
1208 }
1209 //cpl_msg_info(cpl_func,"object_jitters: %p",object_jitters);
1210 //cpl_size num_tot_images = object_jitters->size;
1211 cpl_size edge_llx = 50;
1212 cpl_size edge_lly = 50;
1213 cpl_size edge_urx = 50;
1214 cpl_size edge_ury = 50;
1215 cpl_mask* omask = cpl_mask_new(sx, sy);
1216 cpl_binary* pmask = cpl_mask_get_data(omask);
1217
1218 /* flag left side Y edge */
1219 for(cpl_size j = 0; j < edge_lly; j++) {
1220 for(cpl_size i = 0; i < sx; i++) {
1221 pmask[i+j*sx] = CPL_BINARY_1;
1222 }
1223 }
1224 /* flag right side Y edge */
1225 for(cpl_size j = sy - edge_ury; j < sy; j++) {
1226 for(cpl_size i = 0; i < sx; i++) {
1227 pmask[i+j*sx] = CPL_BINARY_1;
1228 }
1229 }
1230
1231 /* flag bottom side X edge */
1232 for(cpl_size j = 0; j < sy; j++) {
1233 for(cpl_size i = 0; i < edge_llx; i++) {
1234 pmask[i+j*sx] = CPL_BINARY_1;
1235 }
1236 }
1237 /* flag upper side X edge */
1238 for(cpl_size j = 0; j < sy; j++) {
1239 for(cpl_size i = sx - edge_urx; i < sx; i++) {
1240 pmask[i+j*sx] = CPL_BINARY_1;
1241 }
1242 }
1243
1244 for (cpl_size k = 0; k < num_tot_images; k++) {
1245
1246 const hdrl_image * himg = hdrl_imagelist_get_const(science_iml, k);
1247 hdrl_value med_img = hdrl_image_get_median(himg);
1248 hdrl_value scale = {0, 0};
1249 scale.data = med_img.data / med_sky.data;
1250 cpl_msg_warning(cpl_func,"Final frame: %lld med_img.data: %g med_sky_0.data: %g scale.data: %g",
1251 k, med_img.data, med_sky.data, scale.data);
1252 hdrl_image * skysub = hdrl_image_duplicate(himg);
1253 hdrl_image * scaled_sky = hdrl_image_duplicate(sky_flat);
1254 hdrl_image_mul_scalar(scaled_sky, scale);
1255 hdrl_image_sub_image(skysub, scaled_sky);
1256 hdrl_value med_skyub = hdrl_image_get_median(skysub);
1257 /* Create output filename */
1258 char outfile[256];
1259 snprintf(outfile, 256, "sci_skysub_2nd_%03lld.fits", k);
1260
1261 /* Load and update header */
1262 cpl_msg_info(cpl_func,"ESO QC SKY_MED: %g",med_img.data);
1263 cpl_msg_info(cpl_func,"ESO QC SKYSUB_MED: %g", med_skyub.data);
1264
1265 /*
1266
1267 cpl_frame * frame = cpl_frameset_get_position(raw_frames, i);
1268 cpl_propertylist* plist = cpl_propertylist_load(cpl_frame_get_filename(frame), 0);
1269 cpl_propertylist* dhead = cpl_propertylist_load(cpl_frame_get_filename(frame), 1);
1270 cpl_propertylist* ehead = cpl_propertylist_load(cpl_frame_get_filename(frame), 2);
1271 cpl_propertylist* qhead = cpl_propertylist_load(cpl_frame_get_filename(frame), 3);
1272 cpl_propertylist_save(plist, outfile, CPL_IO_CREATE);
1273 cpl_image_save(hdrl_image_get_image(skysub), outfile, CPL_TYPE_DOUBLE, dhead, CPL_IO_EXTEND);
1274 cpl_image_save(hdrl_image_get_error(skysub), outfile, CPL_TYPE_DOUBLE, ehead, CPL_IO_EXTEND);
1275 cpl_mask_save(hdrl_image_get_mask(skysub), outfile, qhead, CPL_IO_EXTEND);
1276 cpl_propertylist_delete(dhead);
1277 cpl_propertylist_delete(ehead);
1278 cpl_propertylist_delete(qhead);
1279 */
1280
1281 if(object_jitters->limages[k]->himage) {
1282 hdrl_image_delete(object_jitters->limages[k]->himage);
1283 object_jitters->limages[k]->himage=hdrl_image_duplicate(skysub);
1284 }
1285
1286 if(object_jitters->limages[k]->bkg) {
1287 hdrl_image_delete(object_jitters->limages[k]->bkg);
1288 object_jitters->limages[k]->bkg=hdrl_image_duplicate(scaled_sky);
1289 }
1290
1291 if(object_jitters->limages[k]->confidence) {
1292 cpl_image_delete(object_jitters->limages[k]->confidence);
1293 object_jitters->limages[k]->confidence=cpl_image_duplicate(confidence);
1294 }
1295 int* pconf = cpl_image_get_data_int(object_jitters->limages[k]->confidence);
1296 object_jitters->limages[k]->object_mask = cpl_mask_duplicate(omask);
1297
1298 /* flag left side Y edge */
1299 for(cpl_size j = 0; j < edge_lly; j++) {
1300 for(cpl_size i = 0; i < sx; i++) {
1301 pconf[i+j*sx] = 0;
1302 }
1303 }
1304 /* flag right side Y edge */
1305 for(cpl_size j = sy - edge_ury; j < sy; j++) {
1306 for(cpl_size i = 0; i < sx; i++) {
1307 pconf[i+j*sx] = 0;
1308 }
1309 }
1310
1311 /* flag bottom side X edge */
1312 for(cpl_size j = 0; j < sy; j++) {
1313 for(cpl_size i = 0; i < edge_llx; i++) {
1314 pconf[i+j*sx] = 0;
1315 }
1316 }
1317 /* flag upper side X edge */
1318 for(cpl_size j = 0; j < sy; j++) {
1319 for(cpl_size i = sx - edge_urx; i < sx; i++) {
1320 pconf[i+j*sx] = 0;
1321 }
1322 }
1323
1324 hdrl_image_delete(skysub);
1325 hdrl_image_delete(scaled_sky);
1326 }
1327 cpl_mask_delete(omask);
1328 eris_check_error_code(cpl_func);
1329 return object_jitters;
1330}
1331
1332
1333
1334static located_imagelist *
1335eris_create_final_skysub_products_robert(const cpl_frameset* frameset, const cpl_frameset* raw_frames, const hdrl_image* sky_flat,
1336 cpl_image* confidence, cpl_frameset** products)
1337{
1338
1339 cpl_msg_info(cpl_func, "Saving sky-subtracted frames...");
1340 cpl_size num_frames = cpl_frameset_get_size(raw_frames);
1341
1342 hdrl_value med_sky = hdrl_image_get_median(sky_flat);
1343 located_imagelist* object_jitters = NULL;
1344
1345 cpl_size sx = hdrl_image_get_size_x(sky_flat);
1346 cpl_size sy = hdrl_image_get_size_y(sky_flat);
1347
1348
1349 cpl_mask* omask = cpl_mask_new(sx, sy);
1350 cpl_binary* pmask = cpl_mask_get_data(omask);
1351
1352 located_image * limage = NULL;
1353 cpl_image * copyconf = NULL;
1354 cpl_size collapse_cube = 0;
1355 cpl_frameset_iterator * frameset_iter = cpl_frameset_iterator_new(raw_frames);
1356 cpl_size k = 0;
1357 const char* fname;
1358 cpl_propertylist* plist = NULL;
1359 cpl_propertylist* hdata = NULL;
1360 cpl_propertylist* herrs = NULL;
1361 cpl_propertylist* hqual = NULL;
1362 cpl_size ext_prim = 0;
1363 cpl_size ext_data = 1;
1364 cpl_size ext_errs = 2;
1365 cpl_size ext_qual = 3;
1366 cpl_image* data = NULL;
1367 cpl_image* errs = NULL;
1368 cpl_image* qual = NULL;
1369 hdrl_image* hima = NULL;
1370
1371 /*
1372 if(cpl_frameset_count_tags(frameset, ERIS_NIX_CAL_DET_OBJECT_JITTER_PRO_CATG) > 0) {
1373 object_jitters = enu_limlist_load_from_frameset(frameset,
1374 ERIS_NIX_CAL_DET_OBJECT_JITTER_PRO_CATG, raw_frames);
1375 } else if(cpl_frameset_count_tags(frameset, ERIS_NIX_CAL_DET_STD_JITTER_PRO_CATG) > 0) {
1376 object_jitters = enu_limlist_load_from_frameset(frameset,
1377 ERIS_NIX_CAL_DET_STD_JITTER_PRO_CATG, raw_frames);
1378 }
1379 */
1380
1381
1382
1383 for (cpl_size iframe = 0; iframe < num_frames; iframe++) {
1384 cpl_frame * frame = cpl_frameset_iterator_get(frameset_iter);
1385 fname = cpl_frame_get_filename(frame);
1386 plist = cpl_propertylist_load(fname, ext_prim);
1387 hdata = cpl_propertylist_load(fname, ext_data);
1388 herrs = cpl_propertylist_load(fname, ext_errs);
1389 hqual = cpl_propertylist_load(fname, ext_qual);
1390 cpl_imagelist* data_iml = cpl_imagelist_load(fname, CPL_TYPE_DOUBLE, ext_data);
1391 cpl_imagelist* errs_iml = cpl_imagelist_load(fname, CPL_TYPE_DOUBLE, ext_errs);
1392 cpl_mask* qual_mask = cpl_mask_load(fname, 0, ext_qual);
1393
1394 cpl_msg_info(cpl_func, "processing %s", cpl_frame_get_filename(frame));
1395 hdrl_imagelist* hlist = hdrl_imagelist_new();
1396 cpl_imagelist* data_list = cpl_imagelist_new();
1397 cpl_imagelist* errs_list = cpl_imagelist_new();
1398 cpl_imagelist* qual_list = cpl_imagelist_new();
1399 cpl_size sz = cpl_imagelist_get_size(data_iml);
1400 cpl_msg_info(cpl_func, "number of planes %lld", sz);
1401 char* outfile = NULL;
1402 for(k = 0; k < sz; k++) {
1403 data = cpl_imagelist_get(data_iml, k);
1404 errs = cpl_imagelist_get(errs_iml, k);
1405 cpl_image_reject_from_mask(data, qual_mask);
1406 cpl_image_reject_from_mask(errs, qual_mask);
1407 hima = hdrl_image_create(data, errs);
1408 hdrl_value med_img = hdrl_image_get_median(hima);
1409 hdrl_value scale = {0, 0};
1410 scale.data = med_img.data / med_sky.data;
1411 cpl_msg_warning(cpl_func,"Final frame: %lld plane %lld med_img.data: %g med_sky_0.data: %g scale.data: %g",
1412 iframe, k, med_img.data, med_sky.data, scale.data);
1413 hdrl_image * skysub = hdrl_image_duplicate(hima);
1414 hdrl_image * scaled_sky = hdrl_image_duplicate(sky_flat);
1415 hdrl_image_mul_scalar(scaled_sky, scale);
1416 hdrl_image_sub_image(skysub, scaled_sky);
1417
1418 hdrl_value med_skyub = hdrl_image_get_median(skysub);
1419
1420 /*
1421 if(object_jitters->limages[k]->himage) {
1422 hdrl_image_delete(object_jitters->limages[k]->himage);
1423 object_jitters->limages[k]->himage=hdrl_image_duplicate(skysub);
1424 }
1425
1426 if(object_jitters->limages[k]->bkg) {
1427 hdrl_image_delete(object_jitters->limages[k]->bkg);
1428 object_jitters->limages[k]->bkg=hdrl_image_duplicate(scaled_sky);
1429 }
1430
1431 if(object_jitters->limages[k]->confidence) {
1432 cpl_image_delete(object_jitters->limages[k]->confidence);
1433 object_jitters->limages[k]->confidence=cpl_image_duplicate(confidence);
1434 }
1435 int* pconf = cpl_image_get_data_int(object_jitters->limages[k]->confidence);
1436 object_jitters->limages[k]->object_mask = cpl_mask_duplicate(omask);
1437 */
1438
1439 /* Create output filename */
1440
1441
1442
1443 /* Load and update header */
1444 cpl_msg_info(cpl_func,"ESO QC SKY_MED: %g",med_img.data);
1445 cpl_msg_info(cpl_func,"ESO QC SKYSUB_MED: %g", med_skyub.data);
1446 if(debug) {
1447 if(k == 0) {
1448 outfile = cpl_sprintf("sci_skysub_2nd_%03lld.fits", iframe);
1449 cpl_propertylist_save(plist, outfile, CPL_IO_CREATE);
1450 } else {
1451 cpl_image_save(hdrl_image_get_image(skysub), outfile, CPL_TYPE_DOUBLE, hdata, CPL_IO_EXTEND);
1452 cpl_image_save(hdrl_image_get_error(skysub), outfile, CPL_TYPE_DOUBLE, herrs, CPL_IO_EXTEND);
1453 //cpl_mask_save(hdrl_image_get_mask(skysub), outfile, hqual, CPL_IO_EXTEND);
1454 }
1455 }
1456 hdrl_imagelist_set(hlist, skysub, k);
1457
1458 hdrl_image_delete(scaled_sky);
1459 hdrl_image_delete(hima);
1460
1461 }
1462
1463 cpl_free(outfile);
1464
1465 outfile = cpl_sprintf("sci_pro_skysub_2nd_%03lld.fits", iframe);
1466
1467 for (cpl_size i = 0; i < sz; i++) {
1468 cpl_imagelist_set(data_list, hdrl_image_get_image(hdrl_imagelist_get(hlist, i)), i);
1469 cpl_imagelist_set(errs_list, hdrl_image_get_error(hdrl_imagelist_get(hlist, i)), i);
1470 cpl_image * bpm = cpl_image_new_from_mask(hdrl_image_get_mask(hdrl_imagelist_get(hlist, i)));
1471 cpl_imagelist_set(qual_list, bpm, i);
1472
1473 }
1474
1475 cpl_propertylist_save(plist, outfile, CPL_IO_CREATE);
1476 cpl_imagelist_save(data_list, outfile, CPL_TYPE_DOUBLE, hdata, CPL_IO_EXTEND);
1477 cpl_imagelist_save(errs_list, outfile, CPL_TYPE_DOUBLE, herrs, CPL_IO_EXTEND);
1478 cpl_imagelist_save(qual_list, outfile, CPL_TYPE_INT, herrs, CPL_IO_EXTEND);
1479 cpl_mask_save(qual_mask, outfile, hqual, CPL_IO_EXTEND);
1480 cpl_image_save(hdrl_image_get_image(sky_flat), outfile, CPL_TYPE_DOUBLE, NULL, CPL_IO_EXTEND);
1481 cpl_image_save(confidence, outfile, CPL_TYPE_DOUBLE, NULL, CPL_IO_EXTEND);
1482
1483 cpl_frame * product = cpl_frame_new();
1484 cpl_frame_set_filename(product, outfile);
1485 cpl_frame_set_tag(product, ERIS_NIX_IMG_SUPERSKY_SKYFLAT);
1486 cpl_frame_set_type(product, CPL_FRAME_TYPE_IMAGE);
1487 cpl_frame_set_group(product, CPL_FRAME_GROUP_PRODUCT);
1488 cpl_frame_set_level(product, CPL_FRAME_LEVEL_FINAL);
1489 cpl_frameset_insert(*products, product);
1490
1491
1492 cpl_propertylist_delete(plist);
1493 cpl_propertylist_delete(hdata);
1494 cpl_propertylist_delete(herrs);
1495 cpl_propertylist_delete(hqual);
1496
1497 cpl_imagelist_delete(data_iml);
1498 cpl_imagelist_delete(errs_iml);
1499 cpl_mask_delete(qual_mask);
1500 hdrl_imagelist_delete(hlist);
1501
1502
1503 for (cpl_size i = sz-1; i >= 0; i--) {
1504 cpl_imagelist_unset(data_list, i);
1505 cpl_imagelist_unset(errs_list, i);
1506 //cpl_imagelist_unset(qual_list, i);
1507 }
1508
1509 cpl_imagelist_delete(data_list);
1510 cpl_imagelist_delete(errs_list);
1511 cpl_imagelist_delete(qual_list);
1512 //cpl_msg_info(cpl_func,"about to exit");
1513 //exit(0);
1514
1515 cpl_free(outfile);
1516
1517
1518 }
1519
1520 cpl_frameset_iterator_delete(frameset_iter);
1521 cpl_mask_delete(omask);
1522 eris_check_error_code(cpl_func);
1523
1524 return object_jitters;
1525}
1526
1527
1528
1529/*----------------------------------------------------------------------------*/
1539/*----------------------------------------------------------------------------*/
1540static cpl_error_code
1541eris_load_and_error_cube(const char* fname, const cpl_mask* bad_pixel_map, cpl_size* kkk,
1542 hdrl_imagelist** science_images){
1543
1544 eris_print_rec_status(0);
1545 cpl_imagelist * data_iml = cpl_imagelist_load(fname, CPL_TYPE_DOUBLE, 1);
1546 cpl_imagelist * error_iml = cpl_imagelist_load(fname, CPL_TYPE_DOUBLE, 2);
1547 cpl_mask * dqual_msk = cpl_mask_load(fname, 0, 3);
1548 cpl_imagelist * confm_iml = cpl_imagelist_load(fname, CPL_TYPE_DOUBLE, 4);
1549
1550 double* pconf_img = cpl_image_get_data_double(cpl_imagelist_get(confm_iml,0));
1551
1552 cpl_binary* pbpm = cpl_mask_get_data(dqual_msk);
1553 cpl_size sx = cpl_mask_get_size_x(dqual_msk);
1554 cpl_size sy = cpl_mask_get_size_y(dqual_msk);
1555 cpl_size sz = cpl_imagelist_get_size(data_iml);
1556
1557 cpl_mask* mask_tot = cpl_mask_duplicate(dqual_msk);
1558 if(bad_pixel_map) {
1559 cpl_mask_or(mask_tot, bad_pixel_map);
1560 }
1561 double* pdata = NULL;
1562 cpl_image * data_img = NULL;
1563 cpl_image * error_img = NULL;
1564 cpl_size pixel = 0;
1565 cpl_size j_sx = 0;
1566 for (cpl_size k = 0; k < sz; k++) {
1567 data_img = cpl_imagelist_get(data_iml,k);
1568 error_img = cpl_imagelist_get(error_iml,k);
1569 pdata = cpl_image_get_data(cpl_imagelist_get(data_iml,k));
1570 for (cpl_size j = 0; j < sy; j++) {
1571
1572 j_sx = j * sx;
1573 for (cpl_size i = 0; i < sx; i++) {
1574
1575 pixel = i + j_sx;
1576 if((pconf_img[pixel] == 0) || !isfinite(pdata[pixel])) {
1577 pbpm[pixel] = CPL_BINARY_1;
1578 }
1579 }
1580
1581 }
1582
1583 cpl_image_reject_from_mask(data_img, mask_tot);
1584 cpl_image_reject_from_mask(error_img, mask_tot);
1585 hdrl_image * himg = hdrl_image_create(data_img, error_img);
1586
1587 hdrl_imagelist_set(*science_images, himg, *kkk);
1588 (*kkk)++;
1589 }
1590
1591 cpl_mask_delete(mask_tot);
1592 cpl_mask_delete(dqual_msk);
1593 cpl_imagelist_delete(data_iml);
1594 cpl_imagelist_delete(error_iml);
1595 cpl_imagelist_delete(confm_iml);
1596
1597 eris_check_error_code(cpl_func);
1598 return cpl_error_get_code();
1599
1600}
1601
1602
1603
1604/*----------------------------------------------------------------------------*/
1614/*----------------------------------------------------------------------------*/
1615static cpl_error_code
1616eris_load_and_error_cube_slice(const char* fname, const cpl_mask* bad_pixel_map, const cpl_size k,
1617 const cpl_size kkk, hdrl_imagelist** science_images){
1618
1619
1620 cpl_imagelist * data_iml = cpl_imagelist_load(fname, CPL_TYPE_DOUBLE, 1);
1621 cpl_imagelist * error_iml = cpl_imagelist_load(fname, CPL_TYPE_DOUBLE, 2);
1622 cpl_mask * dqual_msk = cpl_mask_load(fname, 0, 3);
1623 cpl_imagelist * confm_iml = cpl_imagelist_load(fname, CPL_TYPE_DOUBLE, 4);
1624
1625 double* pconf_img = cpl_image_get_data_double(cpl_imagelist_get(confm_iml,0));
1626
1627 cpl_binary* pbpm = cpl_mask_get_data(dqual_msk);
1628 cpl_size sx = cpl_mask_get_size_x(dqual_msk);
1629 cpl_size sy = cpl_mask_get_size_y(dqual_msk);
1630
1631 cpl_mask* mask_tot = cpl_mask_duplicate(dqual_msk);
1632 if(bad_pixel_map) {
1633 cpl_mask_or(mask_tot, bad_pixel_map);
1634 }
1635 double* pdata = NULL;
1636 cpl_image * data_img = NULL;
1637 cpl_image * error_img = NULL;
1638 cpl_size pixel = 0;
1639 cpl_size j_sx = 0;
1640
1641 data_img = cpl_imagelist_get(data_iml,k);
1642 error_img = cpl_imagelist_get(error_iml,k);
1643 pdata = cpl_image_get_data(cpl_imagelist_get(data_iml,k));
1644
1645 for (cpl_size j = 0; j < sy; j++) {
1646
1647 j_sx = j * sx;
1648 for (cpl_size i = 0; i < sx; i++) {
1649
1650 pixel = i + j_sx;
1651 if((pconf_img[pixel] == 0) || !isfinite(pdata[pixel])) {
1652 pbpm[pixel] = CPL_BINARY_1;
1653 }
1654 }
1655
1656 }
1657
1658 cpl_image_reject_from_mask(data_img, mask_tot);
1659 cpl_image_reject_from_mask(error_img, mask_tot);
1660 hdrl_image * himg = hdrl_image_create(data_img, error_img);
1661 hdrl_imagelist_set(*science_images, himg, kkk);
1662
1663 cpl_mask_delete(mask_tot);
1664 cpl_mask_delete(dqual_msk);
1665 cpl_imagelist_delete(data_iml);
1666 cpl_imagelist_delete(error_iml);
1667 cpl_imagelist_delete(confm_iml);
1668
1669 eris_check_error_code(cpl_func);
1670 return cpl_error_get_code();
1671
1672}
1673
1674
1675
1676
1677
1678static cpl_error_code
1679eris_load_and_error_cube_robert(const char* fname, const cpl_mask* bad_pixel_map,
1680 const char* method, const double kappa, const int niter, const cpl_size kkk,
1681 hdrl_imagelist** science_images){
1682
1683
1684 cpl_imagelist * data_iml = cpl_imagelist_load(fname, CPL_TYPE_DOUBLE, 1);
1685 cpl_imagelist * error_iml = cpl_imagelist_load(fname, CPL_TYPE_DOUBLE, 2);
1686 cpl_mask * dqual_msk = cpl_mask_load(fname, 0, 3);
1687 cpl_imagelist * confm_iml = cpl_imagelist_load(fname, CPL_TYPE_DOUBLE, 4);
1688
1689 double* pconf_img = cpl_image_get_data_double(cpl_imagelist_get(confm_iml,0));
1690
1691 cpl_binary* pbpm = cpl_mask_get_data(dqual_msk);
1692
1693 cpl_size sx = cpl_mask_get_size_x(dqual_msk);
1694 cpl_size sy = cpl_mask_get_size_y(dqual_msk);
1695 cpl_size sz = cpl_imagelist_get_size(data_iml);
1696
1697 cpl_mask* mask_tot = cpl_mask_duplicate(dqual_msk);
1698 if(bad_pixel_map) {
1699 cpl_mask_or(mask_tot, bad_pixel_map);
1700 }
1701 double* pdata = NULL;
1702 cpl_image * data_img = NULL;
1703 cpl_image * error_img = NULL;
1704 cpl_size pixel = 0;
1705 cpl_size j_sx = 0;
1706 eris_print_rec_status(10);
1707 hdrl_imagelist* storage = hdrl_imagelist_new();
1708 for (cpl_size k = 0; k < sz; k++) {
1709
1710 eris_print_rec_status(11);
1711 data_img = cpl_imagelist_get(data_iml,k);
1712 eris_print_rec_status(12);
1713 error_img = cpl_imagelist_get(error_iml,k);
1714 eris_print_rec_status(13);
1715 pdata = cpl_image_get_data(cpl_imagelist_get(data_iml,k));
1716 eris_print_rec_status(14);
1717
1718 for (cpl_size j = 0; j < sy; j++) {
1719
1720 j_sx = j * sx;
1721 for (cpl_size i = 0; i < sx; i++) {
1722
1723 pixel = i + j_sx;
1724 if((pconf_img[pixel] == 0) || !isfinite(pdata[pixel])) {
1725 pbpm[pixel] = CPL_BINARY_1;
1726 }
1727 }
1728
1729 }
1730
1731 cpl_image_reject_from_mask(data_img, mask_tot);
1732 cpl_image_reject_from_mask(error_img, mask_tot);
1733 hdrl_image * himg = hdrl_image_create(data_img, error_img);
1734 hdrl_imagelist_set(storage, himg, k);
1735
1736 }
1737 cpl_mask_delete(mask_tot);
1738 cpl_mask_delete(dqual_msk);
1739 const hdrl_parameter * collapse_par;
1740
1741 if (strcmp(method, "median") == 0) {
1743 } else {
1744 cpl_msg_info(cpl_func,"kappa: %g niter: %d",kappa, niter);
1745 collapse_par = hdrl_collapse_sigclip_parameter_create(kappa, kappa, niter);
1746 }
1747 hdrl_image* image_mean;
1748 cpl_image* conf_map_mean;
1749 hdrl_imagelist_collapse_mean(storage,&image_mean,&conf_map_mean);
1750 hdrl_image* collapsed = NULL;
1751 cpl_image* contribution_map;
1752 hdrl_imagelist_collapse(storage,collapse_par, &collapsed, &contribution_map);
1753
1754 hdrl_imagelist_set(*science_images, image_mean, kkk);
1755 const char* filename = cpl_sprintf("sci_mean_%lld.fits", kkk);
1756 cpl_image_save(hdrl_image_get_image(image_mean), filename, CPL_TYPE_DOUBLE, NULL, CPL_IO_DEFAULT);
1757 free(filename);
1758
1759 filename = cpl_sprintf("sci_mean_clip_%lld.fits", kkk);
1760 cpl_image_save(hdrl_image_get_image(collapsed), filename, CPL_TYPE_DOUBLE, NULL, CPL_IO_DEFAULT);
1761 free(filename);
1762
1763 hdrl_parameter_delete(collapse_par);
1764 cpl_image_delete(conf_map_mean);
1765 cpl_image_delete(contribution_map);
1766 hdrl_image_delete(collapsed);
1767 cpl_imagelist_delete(data_iml);
1768 cpl_imagelist_delete(error_iml);
1769 cpl_imagelist_delete(confm_iml);
1770
1771 hdrl_imagelist_delete(storage);
1772 eris_check_error_code(cpl_func);
1773 return cpl_error_get_code();
1774
1775}
1776
1777/*----------------------------------------------------------------------------*/
1787/*----------------------------------------------------------------------------*/
1788static cpl_error_code
1789eris_load_data_and_error_simple(const char* fname, const cpl_size kk, const cpl_mask* bad_pixel_map,
1790 hdrl_imagelist** science_images)
1791{
1792
1793 cpl_ensure_code(fname, CPL_ERROR_NULL_INPUT);
1794 cpl_image * data_img = cpl_image_load(fname, CPL_TYPE_DOUBLE, 0, 1);
1795 cpl_image * error_img = cpl_image_load(fname, CPL_TYPE_DOUBLE, 0, 2);
1796 cpl_mask* dqual_msk = cpl_mask_load(fname, 0, 3);
1797 cpl_image* confm_img = cpl_image_load(fname, CPL_TYPE_DOUBLE, 0, 4);
1798 double* pdata = cpl_image_get_data(data_img);
1799 double* pconf_img = cpl_image_get_data_double(confm_img);
1800 cpl_binary* pbpm = cpl_mask_get_data(dqual_msk);
1801 cpl_size sx = cpl_image_get_size_x(confm_img);
1802 cpl_size sy = cpl_image_get_size_y(confm_img);
1803 cpl_size j_sx = 0;
1804 cpl_size pixel = 0;
1805 for (cpl_size j = 0; j < sy; j++) {
1806
1807 j_sx = j * sx;
1808 for (cpl_size i = 0; i < sx; i++) {
1809
1810 pixel = i + j_sx;
1811 if((pconf_img[pixel] == 0) || !isfinite(pdata[pixel])) {
1812 pbpm[pixel] = CPL_BINARY_1;
1813 }
1814 }
1815
1816 }
1817
1818 cpl_mask* mask_tot = cpl_mask_duplicate(dqual_msk);
1819 if(bad_pixel_map) {
1820 cpl_mask_or(mask_tot, bad_pixel_map);
1821 }
1822 cpl_image_reject_from_mask(data_img, mask_tot);
1823 cpl_image_reject_from_mask(error_img, mask_tot);
1824
1825 //if(cpl_error_get_code != CPL_ERROR_NONE) {
1826 // exit(0);
1827 //}
1828 hdrl_image * himg = hdrl_image_create(data_img, error_img);
1829 hdrl_imagelist_set(*science_images, himg, kk);
1830
1831 cpl_image_delete(data_img);
1832 cpl_image_delete(error_img);
1833 cpl_mask_delete(mask_tot);
1834 cpl_mask_delete(dqual_msk);
1835 cpl_image_delete(confm_img);
1836
1837 // determine statistics of the science images
1838 cpl_mask* mask = cpl_image_get_bpm(hdrl_image_get_image(himg));
1839 double median=0, mean=0, stdev=0,tmean=0, tstd=0, mad=0, min_val=0, max_val=0;
1840 eris_image_stats(
1841 himg,
1842 2.0, // l_sig: lower sigma clip
1843 2.0, // u_sig: upper sigma clip
1844 mask, // optional mask
1845 0.1, // 10% trim fraction
1846 &median, &mean, &stdev, &tmean, &tstd, &mad, &min_val, &max_val
1847 );
1848
1849 cpl_msg_info(cpl_func,"process file: %s", fname);
1850 cpl_msg_info(cpl_func, "Median=%.3f, Mean=%.3f, StdDev=%.3f "
1851 "Trimmed: Mean=%.3f, StdDev=%.3f, MAD=%.3f "
1852 "Range: [%.3f, %.3f]", median, mean, stdev, tmean, tstd, mad, min_val, max_val);
1853
1854 eris_check_error_code(cpl_func);
1855 return cpl_error_get_code();
1856}
1857
1858/*----------------------------------------------------------------------------*/
1868/*----------------------------------------------------------------------------*/
1869static cpl_error_code
1870eris_load_data_and_crea_error(const char* fname, const cpl_size kk, const cpl_mask* bad_pixel_map,
1871 hdrl_imagelist** science_images)
1872{
1873 cpl_ensure_code(fname, CPL_ERROR_NULL_INPUT);
1874 cpl_image * data_img = NULL;
1875 cpl_image * error_img = NULL;
1876 cpl_msg_info(cpl_func,"sci data with less than 4 extensions");
1877 data_img = cpl_image_load(fname, CPL_TYPE_DOUBLE, 0, 0);
1878 if (!data_img) {
1879 cpl_msg_error(cpl_func, "Failed to load frame %s", fname);
1880 eris_check_error_code(cpl_func);
1881 return cpl_error_get_code();
1882 }
1883 double* pdata = NULL;
1884 cpl_mask * dqual_msk = NULL;
1885 pdata = cpl_image_get_data(data_img);
1886 dqual_msk = cpl_image_get_bpm(data_img);
1887 cpl_size sx = cpl_image_get_size_x(data_img);
1888 cpl_size sy = cpl_image_get_size_y(data_img);
1889 cpl_binary* pbpm = cpl_mask_get_data(dqual_msk);
1890 cpl_size j_sx = 0;
1891 cpl_size pixel = 0;
1892 for (cpl_size j = 0; j < sy; j++) {
1893 j_sx = j * sx;
1894 for (cpl_size i = 0; i < sx; i++) {
1895 pixel = i + j_sx;
1896 if(!isfinite(pdata[pixel])) {
1897 pbpm[pixel] = CPL_BINARY_1;
1898 }
1899 }
1900
1901 }
1902 cpl_mask* mask_tot = cpl_mask_duplicate(dqual_msk);
1903 if(bad_pixel_map) {
1904 cpl_mask_or(mask_tot, bad_pixel_map);
1905 }
1906 // Apply bad pixel mask
1907 cpl_image_reject_from_mask(data_img, mask_tot);
1908
1909 // Create error image (simple Poisson model for now)
1910 error_img = cpl_image_duplicate(data_img);
1911 cpl_image_abs(error_img);
1912 cpl_image_add_scalar(error_img, 1.0); // Add readnoise^2
1913 cpl_image_power(error_img, 0.5);
1914
1915
1916 hdrl_image * himg = hdrl_image_create(data_img, error_img);
1917 hdrl_imagelist_set(*science_images, himg, kk);
1918
1919 cpl_image_delete(data_img);
1920 cpl_image_delete(error_img);
1921 cpl_mask_delete(mask_tot);
1922
1923 eris_check_error_code(cpl_func);
1924 return cpl_error_get_code();
1925}
1926
1927/*----------------------------------------------------------------------------*/
1938/*----------------------------------------------------------------------------*/
1939static cpl_error_code
1940eris_locate_and_mask_sources(hdrl_imagelist* skysub_images, const cpl_image* bpm_image,
1941 const cpl_parameterlist* parlist, const char* prefix, cpl_mask** source_masks)
1942{
1943 cpl_msg_info(cpl_func, "Detecting sources...");
1944 cpl_size num_tot_images = hdrl_imagelist_get_size(skysub_images);
1945 //cpl_msg_info(cpl_func,"num_tot_images: %lld",num_tot_images);
1946 //source_masks = cpl_calloc(num_tot_images, sizeof(cpl_mask *));
1947
1948 hdrl_parameter *par_cat = hdrl_catalogue_parameter_parse_parlist(parlist, prefix);
1949 const hdrl_catalogue_options opt = HDRL_CATALOGUE_SEGMAP;
1951
1952 if(bpm_image) {
1953 if(debug) {
1954 for (int i = 0; i < num_tot_images; i++) {
1955 hdrl_image * himg = hdrl_imagelist_get(skysub_images, i);
1956 cpl_image * data = hdrl_image_get_image(himg);
1957
1958 cpl_mask *bpm = cpl_mask_threshold_image_create(bpm_image, 0, INT_MAX);
1959 cpl_image_reject_from_mask(data, bpm);
1960 cpl_mask_delete(bpm);
1961
1962 source_masks[i] = eris_detect_sources_hdrl_catalogue(data, par_cat);
1963
1964 //char* fname = cpl_sprintf("msk_%d.fits",i);
1965 //cpl_mask_save(source_masks[i], fname, NULL, CPL_IO_DEFAULT);
1966 //cpl_free(fname);
1967
1968 /* Convert BPM image to mask (0=good, >0=bad) */
1969 cpl_mask* bad_pixel_map = cpl_mask_threshold_image_create(bpm_image, 0.5, DBL_MAX);
1970
1971 int nsources = cpl_mask_count(source_masks[i]) - cpl_mask_count(bad_pixel_map);
1972 cpl_msg_info(cpl_func, " Frame %d: detected %d source pixels", i+1, nsources);
1973
1974 cpl_mask_delete(bad_pixel_map);
1975 }
1976 } else {
1977 for (int i = 0; i < num_tot_images; i++) {
1978 hdrl_image * himg = hdrl_imagelist_get(skysub_images, i);
1979 cpl_image * data = hdrl_image_get_image(himg);
1980
1981 cpl_mask *bpm = cpl_mask_threshold_image_create(bpm_image, 0, INT_MAX);
1982 cpl_image_reject_from_mask(data, bpm);
1983 cpl_mask_delete(bpm);
1984
1985 source_masks[i] = eris_detect_sources_hdrl_catalogue(data, par_cat);
1986
1987 /* Convert BPM image to mask (0=good, >0=bad) */
1988 cpl_mask* bad_pixel_map = cpl_mask_threshold_image_create(bpm_image, 0.5, DBL_MAX);
1989
1990 int nsources = cpl_mask_count(source_masks[i]) - cpl_mask_count(bad_pixel_map);
1991 cpl_msg_info(cpl_func, " Frame %d: detected %d source pixels", i+1, nsources);
1992
1993 cpl_mask_delete(bad_pixel_map);
1994 }
1995 }
1996 } else {
1997 if(debug) {
1998 for (int i = 0; i < num_tot_images; i++) {
1999 hdrl_image * himg = hdrl_imagelist_get(skysub_images, i);
2000 cpl_image * data = hdrl_image_get_image(himg);
2001 char* fname = cpl_sprintf("check_data_%d.fits",i);
2002 cpl_image_save(data, fname, CPL_TYPE_DOUBLE, NULL, CPL_IO_DEFAULT);
2003 cpl_image_save(hdrl_image_get_error(himg), fname, CPL_TYPE_DOUBLE, NULL, CPL_IO_EXTEND);
2004 cpl_mask_save(hdrl_image_get_mask(himg), fname, NULL, CPL_IO_EXTEND);
2005 cpl_free(fname);
2006 source_masks[i] = eris_detect_sources_hdrl_catalogue(data, par_cat);
2007 //char* fname = cpl_sprintf("msk_%d.fits",i);
2008 //cpl_mask_save(source_masks[i], fname, NULL, CPL_IO_DEFAULT);
2009 //cpl_free(fname);
2010
2011 int nsources = cpl_mask_count(source_masks[i]);
2012 cpl_msg_info(cpl_func, " Frame %d: detected %d source pixels", i+1, nsources);
2013 }
2014 } else {
2015 for (int i = 0; i < num_tot_images; i++) {
2016 hdrl_image * himg = hdrl_imagelist_get(skysub_images, i);
2017 cpl_image * data = hdrl_image_get_image(himg);
2018
2019 source_masks[i] = eris_detect_sources_hdrl_catalogue(data, par_cat);
2020
2021 int nsources = cpl_mask_count(source_masks[i]);
2022 cpl_msg_info(cpl_func, " Frame %d: detected %d source pixels", i+1, nsources);
2023 }
2024 }
2025 }
2026
2027
2028 if(par_cat) hdrl_parameter_delete(par_cat);
2029 eris_check_error_code(cpl_func);
2030 return cpl_error_get_code();
2031}
2032
2033
2034
2035
2036
2037/*----------------------------------------------------------------------------*/
2048/*----------------------------------------------------------------------------*/
2049static cpl_error_code
2050eris_locate_and_mask_sources2( located_imagelist * jitters, const cpl_image* bpm_image,
2051 const cpl_parameterlist* parlist, const char* prefix, cpl_mask** source_masks)
2052{
2053 cpl_msg_info(cpl_func, "Detecting sources...");
2054 cpl_size num_tot_images = jitters->size;
2055 cpl_msg_info(cpl_func,"num_tot_images3: %lld",num_tot_images);
2056 //source_masks = cpl_calloc(num_tot_images, sizeof(cpl_mask *));
2057
2058 hdrl_parameter *par_cat = hdrl_catalogue_parameter_parse_parlist(parlist, prefix);
2059 const hdrl_catalogue_options opt = HDRL_CATALOGUE_SEGMAP;
2061
2062 if(bpm_image) {
2063 if(debug) {
2064 for (int i = 0; i < num_tot_images; i++) {
2065
2066 cpl_image * data = hdrl_image_get_image(jitters->limages[i]->himage);
2067
2068 cpl_mask *bpm = cpl_mask_threshold_image_create(bpm_image, 0, INT_MAX);
2069 cpl_image_reject_from_mask(data, bpm);
2070 cpl_mask_delete(bpm);
2071
2072 source_masks[i] = eris_detect_sources_hdrl_catalogue(data, par_cat);
2073
2074 //char* fname = cpl_sprintf("msk_%d.fits",i);
2075 //cpl_mask_save(source_masks[i], fname, NULL, CPL_IO_DEFAULT);
2076 //cpl_free(fname);
2077
2078 /* Convert BPM image to mask (0=good, >0=bad) */
2079 cpl_mask* bad_pixel_map = cpl_mask_threshold_image_create(bpm_image, 0.5, DBL_MAX);
2080
2081 int nsources = cpl_mask_count(source_masks[i]) - cpl_mask_count(bad_pixel_map);
2082 cpl_msg_info(cpl_func, " Frame %d: detected %d source pixels", i+1, nsources);
2083
2084 cpl_mask_delete(bad_pixel_map);
2085 }
2086 } else {
2087 for (int i = 0; i < num_tot_images; i++) {
2088
2089 cpl_image * data = hdrl_image_get_image(jitters->limages[i]->himage);
2090
2091 cpl_mask *bpm = cpl_mask_threshold_image_create(bpm_image, 0, INT_MAX);
2092 cpl_image_reject_from_mask(data, bpm);
2093 cpl_mask_delete(bpm);
2094
2095 source_masks[i] = eris_detect_sources_hdrl_catalogue(data, par_cat);
2096
2097 /* Convert BPM image to mask (0=good, >0=bad) */
2098 cpl_mask* bad_pixel_map = cpl_mask_threshold_image_create(bpm_image, 0.5, DBL_MAX);
2099
2100 int nsources = cpl_mask_count(source_masks[i]) - cpl_mask_count(bad_pixel_map);
2101 cpl_msg_info(cpl_func, " Frame %d: detected %d source pixels", i+1, nsources);
2102
2103 cpl_mask_delete(bad_pixel_map);
2104 }
2105 }
2106 } else {
2107 if(debug) {
2108 for (int i = 0; i < num_tot_images; i++) {
2109
2110 cpl_image * data = hdrl_image_get_image(jitters->limages[i]->himage);
2111
2112 source_masks[i] = eris_detect_sources_hdrl_catalogue(data, par_cat);
2113 //char* fname = cpl_sprintf("msk_%d.fits",i);
2114 //cpl_mask_save(source_masks[i], fname, NULL, CPL_IO_DEFAULT);
2115 //cpl_free(fname);
2116
2117 int nsources = cpl_mask_count(source_masks[i]);
2118 cpl_msg_info(cpl_func, " Frame %d: detected %d source pixels", i+1, nsources);
2119 }
2120 } else {
2121 for (int i = 0; i < num_tot_images; i++) {
2122
2123 cpl_image * data = hdrl_image_get_image(jitters->limages[i]->himage);
2124
2125 source_masks[i] = eris_detect_sources_hdrl_catalogue(data, par_cat);
2126
2127 int nsources = cpl_mask_count(source_masks[i]);
2128 cpl_msg_info(cpl_func, " Frame %d: detected %d source pixels", i+1, nsources);
2129 }
2130 }
2131 }
2132
2133
2134 if(par_cat) hdrl_parameter_delete(par_cat);
2135 eris_check_error_code(cpl_func);
2136 return cpl_error_get_code();
2137}
2138
2139
2140/*----------------------------------------------------------------------------*/
2148/*----------------------------------------------------------------------------*/
2149static cpl_error_code
2150eris_first_sky_sub(const hdrl_imagelist* science_images, const hdrl_image* sky_flat_0,
2151 hdrl_imagelist** skysub_images)
2152{
2153
2154 cpl_msg_info(cpl_func, "Subtracting first-pass sky flat...");
2155 *skysub_images = hdrl_imagelist_new();
2156 hdrl_value med_sky_0 = hdrl_image_get_median(sky_flat_0);
2157 char* fname = NULL;
2158 hdrl_value scale = {0, 0};
2159 cpl_size num_tot_images = hdrl_imagelist_get_size(science_images);
2160 for (cpl_size i = 0; i < num_tot_images; i++) {
2161 const hdrl_image * himg = hdrl_imagelist_get_const(science_images, i);
2162 cpl_image_reject_value(hdrl_image_get_image(himg), CPL_VALUE_NOTFINITE);
2163 cpl_image_reject_value(hdrl_image_get_error(himg), CPL_VALUE_NOTFINITE);
2164 //cpl_image_save(hdrl_image_get_image_const(himg),"test.fits", CPL_TYPE_DOUBLE, NULL, CPL_IO_DEFAULT);
2165 //cpl_image_save(hdrl_image_get_error_const(himg),"test.fits", CPL_TYPE_DOUBLE, NULL, CPL_IO_EXTEND);
2166 //cpl_mask_save(hdrl_image_get_mask_const(himg),"test.fits", NULL, CPL_IO_EXTEND);
2167 hdrl_value med_img = hdrl_image_get_median(himg);;
2168
2169 scale.data = med_img.data / med_sky_0.data;
2170 cpl_msg_warning(cpl_func,"1st frame: %lld med_img.data: %g med_sky_0.data: %g scale.data: %g",
2171 i, med_img.data, med_sky_0.data, scale.data);
2172 hdrl_image * skysub = hdrl_image_duplicate(himg);
2173 hdrl_image * scaled_sky = hdrl_image_duplicate(sky_flat_0);
2174 hdrl_image_mul_scalar(scaled_sky, scale);
2175
2176 if(debug) {
2177 fname = cpl_sprintf("sci_%lld.fits",i);
2178 cpl_image_save(hdrl_image_get_image_const(himg),fname, CPL_TYPE_DOUBLE, NULL, CPL_IO_DEFAULT);
2179 cpl_image_save(hdrl_image_get_error_const(himg),fname, CPL_TYPE_DOUBLE, NULL, CPL_IO_EXTEND);
2180 cpl_mask_save(hdrl_image_get_mask_const(himg),fname, NULL, CPL_IO_EXTEND);
2181 cpl_free(fname);
2182 fname = cpl_sprintf("scaled_sky_0_%lld.fits",i);
2183 cpl_image_save(hdrl_image_get_image(scaled_sky),fname, CPL_TYPE_DOUBLE, NULL, CPL_IO_DEFAULT);
2184 cpl_image_save(hdrl_image_get_error(scaled_sky),fname, CPL_TYPE_DOUBLE, NULL, CPL_IO_EXTEND);
2185 cpl_mask_save(hdrl_image_get_mask(scaled_sky),fname, NULL, CPL_IO_EXTEND);
2186 cpl_free(fname);
2187 }
2188
2189 hdrl_image_sub_image(skysub, scaled_sky);
2190
2191 if(debug) {
2192 fname = cpl_sprintf("sci_skysub_1st_%03lld.fits",i);
2193 cpl_image_save(hdrl_image_get_image(skysub), fname, CPL_TYPE_DOUBLE, NULL, CPL_IO_DEFAULT);
2194 cpl_image_save(hdrl_image_get_error(skysub), fname, CPL_TYPE_DOUBLE, NULL, CPL_IO_EXTEND);
2195 cpl_mask_save(hdrl_image_get_mask(skysub), fname, NULL, CPL_IO_EXTEND);
2196 cpl_free(fname);
2197 }
2198 hdrl_imagelist_set(*skysub_images, skysub, i);
2199 hdrl_image_delete(scaled_sky);
2200 //exit(0);
2201 }
2202
2203 eris_check_error_code(cpl_func);
2204 return cpl_error_get_code();
2205}
2206
2207#define HDRL_USE_PRIVATE YES
2208/*----------------------------------------------------------------------------*/
2218/*----------------------------------------------------------------------------*/
2219static hdrl_image * eris_create_sky_flat_first(
2220 hdrl_imagelist * image_list,
2221 cpl_mask * mask,
2222 const char * method,
2223 const double kappa,
2224 const int niter,
2225 const cpl_size plane_id)
2226{
2227 hdrl_image * result = NULL;
2228
2229 hdrl_imagelist * masked_list = NULL;
2230 hdrl_parameter * collapse_par = NULL;
2231 cpl_image* contrib = NULL;
2232
2233 cpl_ensure(image_list, CPL_ERROR_NULL_INPUT, NULL);
2234 cpl_ensure(method, CPL_ERROR_NULL_INPUT, NULL);
2235
2236 int num_images = hdrl_imagelist_get_size(image_list);
2237 cpl_ensure(num_images > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
2238 cpl_msg_warning(cpl_func,"kappa: %g niter: %d", kappa, niter);
2239 /* Apply masks if provided */
2240 if (mask) {
2241 masked_list = hdrl_imagelist_new();
2242
2243 for (int i = 0; i < num_images; i++) {
2244
2245 hdrl_image * himg = hdrl_image_duplicate(
2246 hdrl_imagelist_get_const(image_list, i));
2247
2248 cpl_image * data = hdrl_image_get_image(himg);
2249 cpl_image * errs = hdrl_image_get_error(himg);
2250 cpl_mask * bpm = hdrl_image_get_mask(himg);
2251 cpl_mask_or(mask,bpm);
2252 cpl_image_reject_from_mask(data, mask);
2253 //cpl_image_reject_from_mask(errs, mask);
2254 hdrl_imagelist_set(masked_list, himg, i);
2255
2256 }
2257 } else {
2258 masked_list = image_list;
2259 }
2260
2261 /* Create collapse parameters */
2262 if (strcmp(method, "median") == 0) {
2264 } else {
2265 cpl_msg_info(cpl_func,"kappa: %g niter: %d",kappa, niter);
2266 collapse_par = hdrl_collapse_sigclip_parameter_create(kappa, kappa, niter);
2267 }
2268
2269 /* Collapse image list */
2270 hdrl_imagelist_collapse(image_list, collapse_par, &result, &contrib);
2271 if(debug) {
2272 char* fname;
2273 if(plane_id < 0) {
2274 fname = cpl_sprintf("sky_flat_1st.fits");
2275 } else {
2276 fname = cpl_sprintf("sky_flat_1st_plane_%3.3lld.fits",plane_id);
2277 }
2278 cpl_image_save(hdrl_image_get_image(result), fname,CPL_TYPE_DOUBLE,NULL,CPL_IO_DEFAULT);
2279 cpl_image_save(hdrl_image_get_error(result), fname,CPL_TYPE_DOUBLE,NULL,CPL_IO_EXTEND);
2280 cpl_mask_save(hdrl_image_get_mask(result), fname,NULL,CPL_IO_EXTEND);
2281 cpl_free(fname);
2282 }
2283
2284 /* Cleanup */
2285 cpl_image_delete(contrib);
2286 hdrl_parameter_delete(collapse_par);
2287 if (mask && masked_list) hdrl_imagelist_delete(masked_list);
2288 eris_check_error_code(cpl_func);
2289 return result;
2290}
2291
2292
2293
2294static hdrl_image * eris_create_sky_flat_first2(
2295 located_imagelist * jitters,
2296 located_imagelist ** sky_1st,
2297 cpl_mask * mask,
2298 const char * method,
2299 const double kappa,
2300 const int niter,
2301 const cpl_size plane_id)
2302{
2303 hdrl_image * result = NULL;
2304
2305 //located_imagelist* masked_list = NULL;
2306 hdrl_imagelist * masked_list = NULL;
2307 hdrl_parameter * collapse_par = NULL;
2308 cpl_image* contrib;
2309 cpl_ensure(jitters, CPL_ERROR_NULL_INPUT, NULL);
2310 cpl_ensure(method, CPL_ERROR_NULL_INPUT, NULL);
2311 eris_print_rec_status(0);
2312 int num_images = jitters->size;
2313 cpl_ensure(num_images > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
2314 cpl_msg_warning(cpl_func,"kappa: %g niter: %d", kappa, niter);
2315 /* Apply masks if provided */
2316 eris_print_rec_status(1);
2317 masked_list = hdrl_imagelist_new();
2318 *sky_1st = enu_located_imagelist_duplicate(jitters);
2319 if (mask) {
2320 eris_print_rec_status(2);
2321 for (int i = 0; i < num_images; i++) {
2322 hdrl_image * himg = hdrl_image_duplicate(jitters->limages[i]->himage);
2323 cpl_image * data = hdrl_image_get_image(himg);
2324 cpl_image * errs = hdrl_image_get_error(himg);
2325 cpl_mask * bpm = hdrl_image_get_mask(himg);
2326 cpl_mask_or(mask,bpm);
2327 cpl_image_reject_from_mask(data, mask);
2328 cpl_image_reject_from_mask(errs, mask);
2329 hdrl_imagelist_set(masked_list, himg, i);
2330 }
2331 eris_print_rec_status(3);
2332 } else {
2333 eris_print_rec_status(4);
2334 for (int i = 0; i < num_images; i++) {
2335 //hdrl_image * himg = hdrl_image_duplicate(jitters->limages[i]->himage);
2336 //cpl_image * data = hdrl_image_get_image(himg);
2337 hdrl_imagelist_set(masked_list, jitters->limages[i]->himage, i);
2338 //cpl_image_save(hdrl_image_get_image(himg), "test.fits", CPL_TYPE_DOUBLE, NULL, CPL_IO_DEFAULT);
2339 //cpl_image_save(hdrl_image_get_error(himg), "test.fits", CPL_TYPE_DOUBLE, NULL, CPL_IO_EXTEND);
2340 //cpl_mask_save(hdrl_image_get_mask(himg), "test.fits", NULL, CPL_IO_EXTEND);
2341 //exit(0);
2342 }
2343 eris_print_rec_status(5);
2344 }
2345
2346 /* Create collapse parameters */
2347 if (strcmp(method, "median") == 0) {
2349 } else {
2350 cpl_msg_info(cpl_func,"kappa: %g niter: %d",kappa, niter);
2351 collapse_par = hdrl_collapse_sigclip_parameter_create(kappa, kappa, niter);
2352 }
2353 eris_print_rec_status(6);
2354 /* Collapse image list */
2355 hdrl_imagelist_collapse(masked_list, collapse_par, &result, &contrib);
2356 for (int i = 0; i < num_images; i++) {
2357 hdrl_image_sub_image((*sky_1st)->limages[i]->himage,result);
2358 (*sky_1st)->limages[i]->bkg = hdrl_image_duplicate(result);
2359 (*sky_1st)->limages[i]->confidence = cpl_image_duplicate(contrib);
2360 }
2361 eris_print_rec_status(7);
2362 if(debug) {
2363 char* fname;
2364 if(plane_id < 0) {
2365 fname = cpl_sprintf("sky_flat_1st.fits");
2366 } else {
2367 fname = cpl_sprintf("sky_flat_1st_plane_%3.3lld.fits",plane_id);
2368 }
2369 cpl_image_save(hdrl_image_get_image(result), fname,CPL_TYPE_DOUBLE,NULL,CPL_IO_DEFAULT);
2370 cpl_image_save(hdrl_image_get_error(result), fname,CPL_TYPE_DOUBLE,NULL,CPL_IO_EXTEND);
2371 cpl_mask_save(hdrl_image_get_mask(result), fname,NULL,CPL_IO_EXTEND);
2372 cpl_free(fname);
2373 }
2374 eris_print_rec_status(8);
2375 /* Cleanup */
2376 cpl_image_delete(contrib);
2377 hdrl_parameter_delete(collapse_par);
2378 if (masked_list) {
2379 for (int i = num_images-1; i >= 0; i--) {
2380 hdrl_imagelist_unset(masked_list, i);
2381 }
2382 hdrl_imagelist_delete(masked_list);
2383 }
2384 eris_check_error_code(cpl_func);
2385 return result;
2386}
2387
2388
2398/*----------------------------------------------------------------------------*/
2399static hdrl_image * eris_create_sky_flat_final(
2400 hdrl_imagelist * image_list,
2401 cpl_mask ** masks,
2402 const char * method,
2403 const double kappa,
2404 const int niter,
2405 const cpl_size plane_id,
2406 cpl_image** contrib)
2407{
2408 hdrl_image * result = NULL;
2409 hdrl_imagelist * masked_list = NULL;
2410 hdrl_parameter * collapse_par = NULL;
2411
2412
2413 cpl_ensure(image_list, CPL_ERROR_NULL_INPUT, NULL);
2414 cpl_ensure(method, CPL_ERROR_NULL_INPUT, NULL);
2415
2416 int num_images = hdrl_imagelist_get_size(image_list);
2417 cpl_ensure(num_images > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
2418
2419 /* Apply masks if provided */
2420 if (masks) {
2421 masked_list = hdrl_imagelist_new();
2422
2423 for (int i = 0; i < num_images; i++) {
2424
2425 hdrl_image * himg = hdrl_image_duplicate(
2426 hdrl_imagelist_get_const(image_list, i));
2427
2428 cpl_image * data = hdrl_image_get_image(himg);
2429 cpl_image * errs = hdrl_image_get_image(himg);
2430 cpl_mask * msk = hdrl_image_get_mask(himg);
2431
2432 //AMO: TODO THE FOLLOWING LINE IS SOURCE OF TROUBLES: was
2433 //cpl_image_reject_from_mask(data, masks[i]);
2434 //Now IS (to be sure flag also existing bad pixels)
2435
2436 cpl_mask_or(msk,masks[i]);
2437 cpl_image_reject_from_mask(data, msk);
2438 cpl_image_reject_from_mask(errs, msk);
2439 if(debug) {
2440 char* fname = NULL;
2441 if(plane_id < 0) {
2442 fname = cpl_sprintf("mask_%d.fits",i);
2443 } else {
2444 fname = cpl_sprintf("mask_%d_%3.3lld.fits",i,plane_id);
2445 }
2446 //cpl_msg_warning(cpl_func,"mask[%d]: %p",i, masks[i]);
2447 cpl_mask_save(masks[i], fname, NULL, CPL_IO_DEFAULT);
2448 cpl_free(fname);
2449 }
2450
2451 hdrl_imagelist_set(masked_list, himg, i);
2452
2453 }
2454 } else {
2455 masked_list = image_list;
2456 }
2457
2458 /* Create collapse parameters */
2459 if (strcmp(method, "median") == 0) {
2461 } else {
2462 collapse_par = hdrl_collapse_sigclip_parameter_create(kappa, kappa, niter);
2463 }
2464
2465
2466 /* Collapse image list */
2467 hdrl_imagelist_collapse(masked_list, collapse_par,&result, contrib);
2468 if(debug) {
2469 char* fname;
2470 if(plane_id < 0) {
2471 fname = cpl_sprintf("sky_flat_2nd.fits");
2472 } else {
2473 fname = cpl_sprintf("sky_flat_2nd_plane_%3.3lld.fits",plane_id);
2474 }
2475 cpl_image_save(hdrl_image_get_image(result), fname,CPL_TYPE_DOUBLE,NULL,CPL_IO_DEFAULT);
2476 cpl_image_save(hdrl_image_get_error(result), fname,CPL_TYPE_DOUBLE,NULL,CPL_IO_EXTEND);
2477 cpl_mask_save(hdrl_image_get_mask(result), fname,NULL,CPL_IO_EXTEND);
2478 cpl_free(fname);
2479 }
2480
2481 /* Cleanup */
2482 hdrl_parameter_delete(collapse_par);
2483
2484 if (masks && masked_list) {
2485 hdrl_imagelist_delete(masked_list);
2486 }
2487
2488 eris_check_error_code(cpl_func);
2489 return result;
2490}
2491
2492
2493static hdrl_image * eris_create_sky_flat_final2(
2494 located_imagelist * jitters,
2495 located_imagelist ** sky_2nd,
2496 cpl_mask ** masks,
2497 const char * method,
2498 const double kappa,
2499 const int niter,
2500 const cpl_size plane_id)
2501{
2502 hdrl_image * result = NULL;
2503 hdrl_imagelist * masked_list = NULL;
2504 hdrl_parameter * collapse_par = NULL;
2505 cpl_image* contrib = NULL;
2506
2507 cpl_ensure(jitters, CPL_ERROR_NULL_INPUT, NULL);
2508 cpl_ensure(method, CPL_ERROR_NULL_INPUT, NULL);
2509
2510 cpl_size num_images = jitters->size;
2511 cpl_ensure(num_images > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
2512 cpl_msg_info(cpl_func,"num_images: %lld",num_images);
2513 /* Apply masks if provided */
2514 cpl_msg_warning(cpl_func,"step1");
2515 if (masks) {
2516 masked_list = hdrl_imagelist_new();
2517
2518 for (int i = 0; i < num_images; i++) {
2519
2520 cpl_image * data = hdrl_image_get_image(jitters->limages[i]->himage);
2521
2522 //AMO: TODO THE FOLLOWING LINE IS SOURCE OF TROUBLES
2523 cpl_image_reject_from_mask(data, masks[i]);
2524 if(debug) {
2525 char* fname = NULL;
2526 if(plane_id < 0) {
2527 fname = cpl_sprintf("mask_%d.fits",i);
2528 } else {
2529 fname = cpl_sprintf("mask_%d_%3.3lld.fits",i,plane_id);
2530 }
2531 //cpl_msg_warning(cpl_func,"mask[%d]: %p",i, masks[i]);
2532 cpl_mask_save(masks[i], fname, NULL, CPL_IO_DEFAULT);
2533 cpl_free(fname);
2534 }
2535
2536 hdrl_imagelist_set(masked_list, jitters->limages[i]->himage, i);
2537
2538 }
2539 } else {
2540 for (int i = 0; i < num_images; i++) {
2541 hdrl_imagelist_set(masked_list, jitters->limages[i]->himage, i);
2542 }
2543 }
2544 cpl_msg_warning(cpl_func,"step2");
2545
2546 /* Create collapse parameters */
2547 if (strcmp(method, "median") == 0) {
2549 } else {
2550 collapse_par = hdrl_collapse_sigclip_parameter_create(kappa, kappa, niter);
2551 }
2552 cpl_msg_warning(cpl_func,"step3");
2553
2554 /* Collapse image list */
2555 hdrl_imagelist_collapse(masked_list, collapse_par,&result, &contrib);
2556 if(debug) {
2557 char* fname;
2558 if(plane_id < 0) {
2559 fname = cpl_sprintf("sky_flat_2nd.fits");
2560 } else {
2561 fname = cpl_sprintf("sky_flat_2nd_plane_%3.3lld.fits",plane_id);
2562 }
2563 cpl_image_save(hdrl_image_get_image(result), fname,CPL_TYPE_DOUBLE,NULL,CPL_IO_DEFAULT);
2564 cpl_image_save(hdrl_image_get_error(result), fname,CPL_TYPE_DOUBLE,NULL,CPL_IO_EXTEND);
2565 cpl_mask_save(hdrl_image_get_mask(result), fname,NULL,CPL_IO_EXTEND);
2566 cpl_free(fname);
2567 }
2568 cpl_msg_warning(cpl_func,"step4");
2569
2570 for (int i = 0; i < num_images; i++) {
2571 hdrl_image_sub_image((*sky_2nd)->limages[i]->himage,result);
2572 (*sky_2nd)->limages[i]->confidence = contrib;
2573 (*sky_2nd)->limages[i]->bkg = result;
2574 }
2575 cpl_msg_warning(cpl_func,"step5");
2576
2577 /* Cleanup */
2578 if (contrib != NULL) {
2579 cpl_image_delete(contrib);
2580 }
2581
2582 hdrl_parameter_delete(collapse_par);
2583
2584 if (masks && masked_list) {
2585 for (int i = num_images-1; i >= 0; i--) {
2586 hdrl_imagelist_unset(masked_list,i);
2587 }
2588 hdrl_imagelist_delete(masked_list);
2589 }
2590
2591 eris_check_error_code(cpl_func);
2592 return result;
2593}
2594
2595static cpl_mask * eris_detect_sources_hdrl_catalogue(
2596 const cpl_image * image,
2597 hdrl_parameter* p)
2598{
2599
2600 cpl_mask * mask_res = NULL;
2601
2602 hdrl_catalogue_result *res = hdrl_catalogue_compute(image, NULL, NULL, p);
2603 /*
2604 if(debug) {
2605 char* fname = cpl_sprintf("cat_%d.fits",0);
2606 cpl_table_save (res->catalogue, NULL, NULL, out_fname, CPL_IO_CREATE);
2607 cpl_free(fname);
2608 char* fname = cpl_sprintf("smap_%d.fits",0);
2609 cpl_image_save(res->segmentation_map, fname, CPL_TYPE_DOUBLE, NULL, CPL_IO_DEFAULT);
2610 cpl_free(fname);
2611 char* fname = cpl_sprintf("bkg_%d.fits",0);
2612 cpl_image_save(res->background, fname, CPL_TYPE_DOUBLE, NULL, CPL_IO_DEFAULT);
2613 cpl_free(fname);
2614 }
2615 */
2616 double max = cpl_image_get_max(res->segmentation_map);
2617 mask_res = cpl_mask_threshold_image_create(res->segmentation_map, 0.9, max + 0.1);
2619 eris_check_error_code(cpl_func);
2620 return mask_res;
2621}
2622
2623/*----------------------------------------------------------------------------*/
2630/*----------------------------------------------------------------------------*/
2631static cpl_mask * eris_dilate_mask(
2632 const cpl_mask * input_mask,
2633 int radius)
2634{
2635 cpl_mask * result = NULL;
2636 cpl_mask * temp = NULL;
2637 cpl_size nx, ny;
2638
2639 cpl_ensure(input_mask, CPL_ERROR_NULL_INPUT, NULL);
2640 cpl_ensure(radius > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
2641
2642 nx = cpl_mask_get_size_x(input_mask);
2643 ny = cpl_mask_get_size_y(input_mask);
2644
2645 result = cpl_mask_duplicate(input_mask);
2646
2647 /* Simple box dilation */
2648 for (int iter = 0; iter < radius; iter++) {
2649 temp = cpl_mask_duplicate(result);
2650
2651 for (cpl_size j = 2; j < ny - 1; j++) {
2652 for (cpl_size i = 2; i < nx - 1; i++) {
2653 if (cpl_mask_get(temp, i, j) ||
2654 cpl_mask_get(temp, i-1, j) ||
2655 cpl_mask_get(temp, i+1, j) ||
2656 cpl_mask_get(temp, i, j-1) ||
2657 cpl_mask_get(temp, i, j+1)) {
2658 cpl_mask_set(result, i, j, CPL_BINARY_1);
2659 }
2660 }
2661 }
2662 cpl_mask_delete(temp);
2663 }
2664 eris_check_error_code(cpl_func);
2665 return result;
2666}
2667
2668/*----------------------------------------------------------------------------*/
2698/*----------------------------------------------------------------------------*/
2699cpl_error_code
2700eris_image_stats(const hdrl_image *in_image,
2701 double l_sig,
2702 double u_sig,
2703 const cpl_mask *mask,
2704 double trim_fraction,
2705 double *median,
2706 double *mean,
2707 double *stdev,
2708 double *tmean,
2709 double *tstd,
2710 double *mad,
2711 double *min_val,
2712 double *max_val)
2713{
2714 cpl_ensure_code(in_image, CPL_ERROR_NULL_INPUT);
2715 cpl_ensure_code(median, CPL_ERROR_NULL_INPUT);
2716 cpl_ensure_code(mean, CPL_ERROR_NULL_INPUT);
2717 cpl_ensure_code(stdev, CPL_ERROR_NULL_INPUT);
2718 cpl_ensure_code(tmean, CPL_ERROR_NULL_INPUT);
2719 cpl_ensure_code(tstd, CPL_ERROR_NULL_INPUT);
2720 cpl_ensure_code(mad, CPL_ERROR_NULL_INPUT);
2721 cpl_ensure_code(min_val, CPL_ERROR_NULL_INPUT);
2722 cpl_ensure_code(max_val, CPL_ERROR_NULL_INPUT);
2723
2724 const cpl_image *img = hdrl_image_get_image_const(in_image);
2725 cpl_image *work_img = NULL;
2726 cpl_vector *data_vec = NULL;
2727
2728 /* Create working copy of image */
2729 work_img = cpl_image_duplicate(img);
2730 if (!work_img) {
2731 return cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_OUTPUT,
2732 "Failed to duplicate image");
2733 }
2734
2735 /* Apply mask if provided */
2736 if (mask != NULL) {
2737 cpl_image_reject_from_mask(work_img, mask);
2738 }
2739
2740 /* Step 1: Compute initial statistics (NaN-aware) */
2741 hdrl_value median0 = hdrl_image_get_median(in_image);
2742 hdrl_value mean0 = hdrl_image_get_mean(in_image);
2743 double stdev0 = hdrl_image_get_stdev(in_image);
2744
2745 *median = median0.data;
2746 *mean = mean0.data;
2747 *stdev = stdev0;
2748
2749 /* Step 2: Convert image to vector for trimmed statistics */
2750 cpl_size nx = cpl_image_get_size_x(work_img);
2751 cpl_size ny = cpl_image_get_size_y(work_img);
2752 cpl_size npix = nx * ny;
2753
2754 /* Count valid (non-rejected) pixels */
2755 const cpl_mask *rej_mask = cpl_image_get_bpm_const(work_img);
2756 cpl_size nvalid = npix;
2757 if (rej_mask) {
2758 nvalid = npix - cpl_mask_count(rej_mask);
2759 }
2760
2761 if (nvalid == 0) {
2762 cpl_image_delete(work_img);
2763 return cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
2764 "No valid pixels in image");
2765 }
2766
2767 /* Extract valid pixel values to vector */
2768 data_vec = cpl_vector_new(nvalid);
2769 cpl_size idx = 0;
2770 for (cpl_size j = 1; j <= ny; j++) {
2771 for (cpl_size i = 1; i <= nx; i++) {
2772 int is_rejected = 0;
2773 double val = cpl_image_get(work_img, i, j, &is_rejected);
2774 if (!is_rejected && isfinite(val)) {
2775 cpl_vector_set(data_vec, idx, val);
2776 idx++;
2777 }
2778 }
2779 }
2780
2781 /* Adjust vector size if some NaNs were found */
2782 if (idx < nvalid) {
2783 cpl_vector *tmp = cpl_vector_extract(data_vec, 0, idx - 1, 1);
2784 cpl_vector_delete(data_vec);
2785 data_vec = tmp;
2786 nvalid = idx;
2787 }
2788
2789 /* Step 3: Compute trimmed mean (trim_fraction from each end) */
2790 cpl_vector_sort(data_vec, CPL_SORT_ASCENDING);
2791 cpl_size trim_count = (cpl_size)(trim_fraction * nvalid);
2792 if (trim_count * 2 >= nvalid) {
2793 trim_count = 0; /* Not enough data for trimming */
2794 }
2795
2796 double tmean0_10 = 0.0;
2797 if (trim_count > 0) {
2798 /* Compute mean of middle (1 - 2*trim_fraction) of data */
2799 cpl_vector *trimmed = cpl_vector_extract(data_vec, trim_count,
2800 nvalid - trim_count - 1, 1);
2801 tmean0_10 = cpl_vector_get_mean(trimmed);
2802 cpl_vector_delete(trimmed);
2803 } else {
2804 tmean0_10 = cpl_vector_get_mean(data_vec);
2805 }
2806
2807 /* Step 4: Define clipping range based on trimmed mean */
2808 double rmin = tmean0_10 - l_sig * (*stdev);
2809 double rmax = tmean0_10 + u_sig * (*stdev);
2810
2811 /* Step 5: Compute trimmed statistics within clipping range */
2812 cpl_size nclipped = 0;
2813 double sum_clipped = 0.0;
2814 double sum_sq_clipped = 0.0;
2815
2816 for (cpl_size i = 0; i < nvalid; i++) {
2817 double val = cpl_vector_get(data_vec, i);
2818 if (val >= rmin && val <= rmax) {
2819 sum_clipped += val;
2820 sum_sq_clipped += val * val;
2821 nclipped++;
2822 }
2823 }
2824
2825 if (nclipped > 0) {
2826 *tmean = sum_clipped / nclipped;
2827
2828 if (nclipped > 1) {
2829 double variance = (sum_sq_clipped - sum_clipped * sum_clipped / nclipped)
2830 / (nclipped - 1);
2831 *tstd = (variance > 0.0) ? sqrt(variance) : 0.0;
2832 } else {
2833 *tstd = 0.0;
2834 }
2835 } else {
2836 /* No data within clipping range - use unclipped statistics */
2837 *tmean = *mean;
2838 *tstd = *stdev;
2839 }
2840
2841 /* Step 6: Compute MAD (Median Absolute Deviation) */
2842 /* MAD = median(|data - median|) */
2843 cpl_vector *abs_dev = cpl_vector_duplicate(data_vec);
2844 cpl_vector_subtract_scalar(abs_dev, *median);
2845 for (cpl_size i = 0; i < nvalid; i++) {
2846 double val = cpl_vector_get(abs_dev, i);
2847 cpl_vector_set(abs_dev, i, fabs(val));
2848 }
2849 cpl_vector_sort(abs_dev, CPL_SORT_ASCENDING);
2850
2851 /*AMO Check MAD computation versus HDRL
2852 double sigma;
2853 cpl_image* ima = cpl_image_wrap_double(nvalid, 1, cpl_vector_get_data(abs_dev));
2854
2855 *mad=cpl_image_get_mad_window(ima, 1, 1, nvalid, 1, &sigma);
2856
2857 cpl_msg_warning(cpl_func,"HDRL mad: %g", *mad);
2858 */
2859
2860 *mad = cpl_vector_get_median(abs_dev);
2861 //cpl_msg_warning(cpl_func,"OTHER mad: %g", *mad);
2862 cpl_vector_delete(abs_dev);
2863
2864 /* Step 7: Get min and max */
2865 *min_val = cpl_vector_get_min(data_vec);
2866 *max_val = cpl_vector_get_max(data_vec);
2867
2868 /* Cleanup */
2869 cpl_vector_delete(data_vec);
2870 cpl_image_delete(work_img);
2871 eris_check_error_code(cpl_func);
2872 return CPL_ERROR_NONE;
2873}
cpl_error_code enu_dfs_save_limage(cpl_frameset *allframes, const cpl_parameterlist *parlist, const cpl_frameset *provenance, const cpl_boolean prov_raw, const located_image *limage, const char *recipe, const cpl_frame *inherit, cpl_propertylist *applist, const char *pipe_id, const char *filename)
Save a located image structure to a MEF.
Definition: eris_nix_dfs.c:909
void enu_located_imagelist_delete(located_imagelist *limlist)
Delete a located_imagelist and its contents.
cpl_error_code enu_located_imagelist_insert(located_imagelist *limlist, located_image *limage, cpl_size position)
Insert a located_image at a specified point in a located_imagelist.
located_image * enu_load_limage_from_frame(const cpl_frame *frame, cpl_image **pcopyconf, const cpl_boolean collapse_cube)
Load components of a located_image from a frame.
located_imagelist * enu_located_imagelist_duplicate(const located_imagelist *limlist)
Make a deep copy of a located_imagelist and its contents.
cpl_error_code enu_get_ra_dec(const cpl_wcs *wcs, double *ra, double *dec)
Get RA and Dec at centre of image with given wcs.
cpl_error_code enu_debug_limlist_save(const int debug, const located_imagelist *limlist, const char *nameroot, const char *recipename, cpl_frameset *frameset, const cpl_parameterlist *parlist, const cpl_frameset *used)
Save a list of intermediate image results for use in debugging.
located_imagelist * enu_limlist_load_from_frameset(cpl_frameset *frameset, const char *tag, cpl_frameset *used)
Load tagged data from a frameset into a located_imagelist.
char * enu_repreface(const char *filename, const char *preface)
Preface a raw filename with a string.
located_imagelist * enu_located_imagelist_new(cpl_size size)
Return a pointer to a new located_imagelist.
const char * eris_pfits_get_frame_format(const cpl_propertylist *plist)
find out the frame format (DET FRAM FORMAT) value
Definition: eris_pfits.c:889
cpl_error_code eris_check_error_code(const char *func_id)
handle CPL errors
Definition: eris_utils.c:56
cpl_frameset * eris_dfs_extract_frames_with_tag(cpl_frameset *input, const char *rtag)
Extract frames of user given tag.
Definition: eris_utils.c:227
hdrl_catalogue_result * hdrl_catalogue_compute(const cpl_image *image_, const cpl_image *confidence_map, const cpl_wcs *wcs, hdrl_parameter *param_)
build object catalog
void hdrl_catalogue_result_delete(hdrl_catalogue_result *result)
delete hdrl parameter result object
hdrl_parameter * hdrl_catalogue_parameter_parse_parlist(const cpl_parameterlist *parlist, const char *prefix)
Parse parameter list to create input parameters for the catalogue.
cpl_error_code hdrl_catalogue_parameter_set_option(hdrl_parameter *par, hdrl_catalogue_options opt)
set result option of catalogue parameter
hdrl_parameter * hdrl_collapse_sigclip_parameter_create(double kappa_low, double kappa_high, int niter)
create a parameter object for sigclipped mean
hdrl_parameter * hdrl_collapse_median_parameter_create(void)
create a parameter object for median
cpl_error_code hdrl_image_sub_image(hdrl_image *self, const hdrl_image *other)
Subtract two images, store the result in the first image.
hdrl_value hdrl_image_get_median(const hdrl_image *self)
computes the median and associated error of an image.
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
double hdrl_image_get_stdev(const hdrl_image *self)
computes the standard deviation of the data of an image
cpl_mask * hdrl_image_get_mask(hdrl_image *himg)
get cpl bad pixel mask from image
Definition: hdrl_image.c:157
cpl_image * hdrl_image_get_error(hdrl_image *himg)
get error as cpl image
Definition: hdrl_image.c:131
hdrl_value hdrl_image_get_mean(const hdrl_image *self)
computes mean pixel value and associated error of an image.
cpl_size hdrl_image_get_size_y(const hdrl_image *self)
return size of Y dimension of image
Definition: hdrl_image.c:540
const cpl_mask * hdrl_image_get_mask_const(const hdrl_image *himg)
get cpl bad pixel mask from image
Definition: hdrl_image.c:175
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_error_const(const hdrl_image *himg)
get error as cpl image
Definition: hdrl_image.c:144
hdrl_image * hdrl_image_create(const cpl_image *image, const cpl_image *error)
create a new hdrl_image from to existing images by copying them
Definition: hdrl_image.c:295
cpl_image * hdrl_image_get_image(hdrl_image *himg)
get data as cpl image
Definition: hdrl_image.c:105
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_imagelist_collapse_mean(const hdrl_imagelist *himlist, hdrl_image **out, cpl_image **contrib)
Mean collapsing of image list.
hdrl_image * hdrl_imagelist_unset(hdrl_imagelist *himlist, cpl_size pos)
Remove an image from an imagelist.
cpl_error_code hdrl_imagelist_set(hdrl_imagelist *himlist, hdrl_image *himg, cpl_size pos)
Insert an image into an imagelist.
void hdrl_imagelist_delete(hdrl_imagelist *himlist)
Free all memory used by a hdrl_imagelist object including the images.
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_imagelist * hdrl_imagelist_new(void)
Create an empty imagelist.
hdrl_image * hdrl_imagelist_get(const hdrl_imagelist *himlist, cpl_size inum)
Get an image from a list of images.
void hdrl_parameter_delete(hdrl_parameter *obj)
shallow delete of a parameter