ERIS Pipeline Reference Manual 1.9.2
eris_ifu_functions.c
1/* $Id: eris_ifu_utils.c,v 1.12 2013-03-25 11:46:49 cgarcia Exp $
2 *
3 * This file is part of the ERIS Pipeline
4 * Copyright (C) 2002,2003 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 * This program is distributed in the hope that it will be useful,
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 * GNU General Public License for more details.
15 *
16 * You should have received a copy of the GNU General Public License
17 * along with this program; if not, write to the Free Software
18 * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
19 */
20
21/*
22 * $Author: cgarcia $
23 * $Date: 2013-03-25 11:46:49 $
24 * $Revision: 1.12 $
25 * $Name: not supported by cvs2svn $
26 */
27
28#ifdef HAVE_CONFIG_H
29#include <config.h>
30#endif
31
32/*-----------------------------------------------------------------------------
33 Includes
34 -----------------------------------------------------------------------------*/
35
36#include <string.h>
37#include <gsl/gsl_interp.h>
38#include "eris_ifu_dfs.h"
39#include "eris_ifu_debug.h"
40#include "eris_ifu_utils.h"
41#include "eris_ifu_functions.h"
42#include "eris_ifu_distortion_static.h"
43#include "eris_ifu_flat_static.h"
44#include "eris_ifu_error.h"
45#include "eris_utils.h"
46
68/*----------------------------------------------------------------------------*/
77/*----------------------------------------------------------------------------*/
78ifsInstrument eris_ifu_get_instrument_frame(cpl_frame *frame)
79{
80 ifsInstrument instrument = UNSET_INSTRUMENT;
81 const char *filename = NULL;
82 cpl_propertylist *header = NULL;
83
84 TRY
85 {
86 if (frame == NULL) {
87 BRK_WITH_ERROR(CPL_ERROR_NULL_INPUT);
88 }
89
90 filename = cpl_frame_get_filename(frame);
91 header = cpl_propertylist_load(filename, 0);
92 instrument = eris_ifu_get_instrument(header);
94 }
95 CATCH
96 {
97 instrument = UNSET_INSTRUMENT;
98 }
99
101 eris_check_error_code("eris_ifu_get_instrument_frame");
102 return instrument;
103}
104
105/*----------------------------------------------------------------------------*/
115/*----------------------------------------------------------------------------*/
116hdrl_image * eris_ifu_raw_hdrl_image(const cpl_image *cplImage)
117{
118 hdrl_image *image = NULL;
119 cpl_image *noiseImage = NULL;
120
121 cpl_ensure(cplImage, CPL_ERROR_NULL_INPUT, NULL);
122
123 TRY
124 {
125 noiseImage = eris_ifu_calc_noise_map(cplImage, 2.10, 0.);
126 image = hdrl_image_create (cplImage, noiseImage);
127 }
128 CATCH
129 {
130 }
131
132 cpl_image_delete(noiseImage);
133 eris_check_error_code("eris_ifu_raw_hdrl_image");
134 return image;
135}
136
137/*----------------------------------------------------------------------------*/
148/*----------------------------------------------------------------------------*/
149cpl_image* eris_ifu_calc_noise_map(const cpl_image *data,
150 double gain,
151 double readnoise)
152{
153 cpl_image *noiseImg = NULL;
154
155 cpl_ensure(data, CPL_ERROR_NULL_INPUT, NULL);
156 cpl_ensure(gain > 0.0, CPL_ERROR_ILLEGAL_INPUT, NULL);
157 cpl_ensure(readnoise >= 0.0, CPL_ERROR_ILLEGAL_INPUT, NULL);
158
159 TRY
160 {
162 noiseImg = cpl_image_duplicate(data));
163
164 /* calculate initial noise estimate
165 sigma = sqrt(counts * gain + readnoise^2) / gain */
166 cpl_image_abs(noiseImg);
167 cpl_image_multiply_scalar(noiseImg, gain);
168 cpl_image_add_scalar(noiseImg, readnoise * readnoise);
169 cpl_image_power(noiseImg, 0.5);
170 cpl_image_divide_scalar(noiseImg, gain);
171 }
172 CATCH
173 {
174 CATCH_MSG();
175 eris_ifu_free_image(&noiseImg);
176 }
177 eris_check_error_code("eris_ifu_calc_noise_map");
178 return noiseImg;
179}
180
181/*----------------------------------------------------------------------------*/
195/*----------------------------------------------------------------------------*/
196hdrl_imagelist* eris_ifu_load_exposure_frameset(const cpl_frameset *frameset,
197 int exposureCorrectionMode) {
198 hdrl_imagelist *imageList = NULL;
199 hdrl_image *tmpImage = NULL;
200 const cpl_frame *frame = NULL;
201 cpl_size frameCnt;
202 cpl_ensure(frameset, CPL_ERROR_NULL_INPUT, NULL);
203
204 TRY
205 {
206 frameCnt = cpl_frameset_get_size(frameset);
207 imageList = hdrl_imagelist_new();
208 for (cpl_size ix=0; ix<frameCnt; ix++) {
209 frame = cpl_frameset_get_position_const(frameset, ix);
211 frame, exposureCorrectionMode, NULL);
212 hdrl_imagelist_set(imageList, tmpImage,
213 hdrl_imagelist_get_size(imageList));
214 }
215 }
216 CATCH
217 {
218 imageList = NULL;
219 }
220 eris_check_error_code("eris_ifu_load_exposure_frameset");
221 return imageList;
222
223}
224
225/*----------------------------------------------------------------------------*/
236/*----------------------------------------------------------------------------*/
237hdrl_image* eris_ifu_load_exposure_frame(const cpl_frame *frame,
238 int exposureCorrectionMode,
239 cpl_image *dqi) {
240
241 hdrl_image *image = NULL;
242 const char *filename;
243
244 cpl_ensure(frame, CPL_ERROR_NULL_INPUT, NULL);
245 TRY
246 {
247 filename = cpl_frame_get_filename(frame);
249 filename, exposureCorrectionMode, dqi);
250 }
251 CATCH
252 {
253 image = NULL;
254 }
255 eris_check_error_code("eris_ifu_load_exposure_frame");
256 return image;
257
258}
259
260/*----------------------------------------------------------------------------*/
261//TODO: why should one set flagged pixels to 1 in the image? Isn't sufficient
262//to flag them? That is also a different convention from the one of using NaN
273/*----------------------------------------------------------------------------*/
274cpl_error_code eris_ifu_add_badpix_border(cpl_image *data, cpl_boolean add_ones,
275 cpl_image *dqi)
276{
277 cpl_error_code err = CPL_ERROR_NONE;
278 cpl_size x = 0;
279 cpl_size y = 0;
280
281 cpl_ensure_code(data, CPL_ERROR_NULL_INPUT);
282
283 TRY
284 {
285 for (x = 1; x <= ERIS_IFU_DETECTOR_SIZE_X; x++) {
286 /* reject lower border*/
287 for (y = 1;
288 y <= ERIS_IFU_DETECTOR_BP_BORDER;
289 y++)
290 {
291 cpl_image_reject(data, x, y);
292 if (dqi != NULL) {
293 cpl_image_set(dqi, x, y, ERIS_DQI_BP);
294 }
295 if (add_ones) {
296 cpl_image_set(data, x, y, 1.);
297 }
298 }
299
300 /* reject upper border*/
301 for (y = ERIS_IFU_DETECTOR_SIZE_Y-ERIS_IFU_DETECTOR_BP_BORDER+1;
302 y <= ERIS_IFU_DETECTOR_SIZE_Y;
303 y++)
304 {
305 cpl_image_reject(data, x, y);
306 if (dqi != NULL) {
307 cpl_image_set(dqi, x, y, ERIS_DQI_BP);
308 }
309 if (add_ones) {
310 cpl_image_set(data, x, y, 1.);
311 }
312 }
313 }
315
316 for (y = ERIS_IFU_DETECTOR_BP_BORDER+1; y <= ERIS_IFU_DETECTOR_SIZE_Y-ERIS_IFU_DETECTOR_BP_BORDER; y++) {
317 /* reject remaining left border*/
318 for (x = 1; x <= ERIS_IFU_DETECTOR_BP_BORDER; x++) {
319 cpl_image_reject(data, x, y);
320 if (dqi != NULL) {
321 cpl_image_set(dqi, x, y, ERIS_DQI_BP);
322 }
323 if (add_ones) {
324 cpl_image_set(data, x, y, 1.);
325 }
326 }
327
328 /* reject remaining right border*/
329 for (x = ERIS_IFU_DETECTOR_SIZE_X-ERIS_IFU_DETECTOR_BP_BORDER+1;
330 x <= ERIS_IFU_DETECTOR_SIZE_X;
331 x++)
332 {
333 cpl_image_reject(data, x, y);
334 if (dqi != NULL) {
335 cpl_image_set(dqi, x, y, ERIS_DQI_BP);
336 }
337 if (add_ones) {
338 cpl_image_set(data, x, y, 1.);
339 }
340 }
341 }
343 }
344 CATCH
345 {
346 err = cpl_error_get_code();
347 }
348 eris_check_error_code("eris_ifu_add_badpix_border");
349 return err;
350}
351
352/*----------------------------------------------------------------------------*/
380/*----------------------------------------------------------------------------*/
381hdrl_image* eris_ifu_load_exposure_file(const char *filename,
382 int exposureCorrectionMode,
383 cpl_image *dqi)
384{
385 hdrl_image *image = NULL;
386 cpl_image *dataImage = NULL;
387 cpl_image *noiseImage = NULL;
388 cpl_propertylist *header = NULL;
389 ifsInstrument instrument = UNSET_INSTRUMENT;
390 double gain;
391 double ron;
392
393 cpl_ensure(filename, CPL_ERROR_NULL_INPUT, NULL);
394
395 TRY
396 {
397 dataImage = cpl_image_load(filename, CPL_TYPE_UNSPECIFIED, 0, 0);
398 eris_ifu_add_badpix_border(dataImage, CPL_FALSE, dqi);
399 eris_ifu_saturation_detection(dataImage, dqi);
400
401 // ==========
402 // correct data image
403 // ==========
404 if ((exposureCorrectionMode & LINE_EXPOSURE_CORRECTION) != 0) {
406 // cpl_image_save(dataImage,"line_corr_img.fits",
407 // CPL_TYPE_FLOAT,NULL,CPL_IO_CREATE);
408 }
409 if ((exposureCorrectionMode & COLUMN_EXPOSURE_CORRECTION) != 0) {
411 // cpl_image_save(dataImage,"col_corr_img.fits",
412 // CPL_TYPE_FLOAT,NULL,CPL_IO_CREATE);
413 }
414
415 // ==========
416 // fill noise image
417 // ==========
418 // get detector gain and readout noise parameters from FITS header
419 header = cpl_propertylist_load(filename, 0);
420 instrument = eris_ifu_get_instrument(header);
422 if (instrument == SPIFFI) {
423 //FITS header keywords are missing for SPIFFI, set defaults
424 gain = 2.1;
425 ron = 0.0;
426 } else if (instrument == SPIFFIER) {
427 gain = cpl_propertylist_get_double(header, FHDR_DET_CHIP_GAIN);
428 ron = cpl_propertylist_get_double(header, FHDR_DET_CHIP_RON);
430 } else {
431 gain = 0.;
432 ron = 0.;
433 }
434 noiseImage = eris_ifu_calc_noise_map(dataImage, gain, ron);
435 image = hdrl_image_create(dataImage, noiseImage);
436
437 }
438 CATCH
439 {
440 CATCH_MSGS();
442 }
443
445 eris_ifu_free_image(&dataImage);
446 eris_ifu_free_image(&noiseImage);
447 eris_check_error_code("eris_ifu_load_exposure_file");
448 return image;
449}
450
451/*----------------------------------------------------------------------------*/
465/*----------------------------------------------------------------------------*/
466cpl_mask* eris_ifu_detect_crh(hdrl_image* image,
467 int exposureCorrectionMode,
468 hdrl_parameter *laCosmicParams,
469 bool maskImage)
470{
471 cpl_mask *mask = NULL;
472 cpl_ensure(image, CPL_ERROR_NULL_INPUT, NULL);
473 cpl_ensure(laCosmicParams, CPL_ERROR_NULL_INPUT, NULL);
474
475
476 TRY
477 {
478 if ((exposureCorrectionMode & COSMIC_RAY_EXPOSURE_DETECTION) &&
479 (laCosmicParams != NULL)) {
480 mask = hdrl_lacosmic_edgedetect(image, laCosmicParams);
481 cpl_msg_info(__func__,"CRH count %lld",cpl_mask_count(mask));
482 if (maskImage) {
484 hdrl_image_reject_from_mask(image, mask));
485 }
486 // if (loadingCnt == 0) {
487 // cpl_mask_save(mask,"crh.fits", NULL, CPL_IO_CREATE);
488 // } else {
489 // cpl_mask_save(mask,"crh.fits", NULL, CPL_IO_EXTEND);
490 // }
491 // loadingCnt++;
492 }
493 }
494 CATCH
495 {
496 CATCH_MSGS();
497 eris_ifu_free_mask(&mask);
498 }
499 eris_check_error_code("eris_ifu_detect_crh");
500 return mask;
501}
502
503/*----------------------------------------------------------------------------*/
516/*----------------------------------------------------------------------------*/
517/* TODO: this function expects a floating point image, should rather work with double
518 * filling of dqi to be understood: with real data works as implemented but code it's strange
519 */
520
521cpl_error_code eris_ifu_saturation_detection(cpl_image *image, cpl_image *dqi) {
522
523 cpl_ensure_code(image, CPL_ERROR_NULL_INPUT);
524 //TODO: you should not hard-code numbers, rather use a #define
525 const float satPixelValue = -45727232.0f;
526 const float satPixelLevel = satPixelValue + 1.f;
527
528 cpl_size nx = cpl_image_get_size_x(image);
529 cpl_size ny = cpl_image_get_size_y(image);
530 const float *data = cpl_image_get_data_float_const(image);
531 cpl_mask *mask = cpl_image_get_bpm(image);
532 cpl_binary *bp = cpl_mask_get_data(mask);
533 //int* qp = cpl_image_get_data_int(dqi);
534
535 int nSat = 0;
536 for (cpl_size iy=0; iy<ny; iy++) {
537 for (cpl_size ix=0; ix<nx; ix++) {
538 if (data[ix+iy*nx] <= satPixelLevel) {
539 if (dqi != NULL) {
540 cpl_image_add_scalar(dqi, ERIS_DQI_SAT);
541 //qp[ix+iy*nx] = ERIS_DQI_SAT;
542 }
543 bp[ix+iy*nx] = BAD_PIX;
544 nSat++;
545 }
546 }
547 }
548 cpl_image_set_bpm(image, mask);
549 //TODO: you should not hard-code numbers, rather use a #define
550 if (nSat > 50) {
551 cpl_msg_debug(__FUNCTION__, "Found %d saturated pixels", nSat);
552 }
553 eris_check_error_code("eris_ifu_saturation_detection");
554 return CPL_ERROR_NONE;
555}
556
557/*----------------------------------------------------------------------------*/
577/*----------------------------------------------------------------------------*/
578cpl_error_code eris_ifu_exposure_line_correction(cpl_image *image)
579{
580 cpl_error_code retVal = CPL_ERROR_NONE;
581 cpl_size nCols;
582 cpl_size nRows;
583 cpl_type type;
584 double *dData = NULL;
585 float *fData = NULL;
586 cpl_vector *borderVector = NULL;
587 double *borderData;
588 double borderMean;
589 double borderMedian;
590 float fborderMean;
591 float fborderMedian;
592 const int borderSize = 4;
593
594 cpl_ensure_code(image, CPL_ERROR_NULL_INPUT);
595 TRY
596 {
598 borderVector = cpl_vector_new(borderSize*2));
599 borderData = cpl_vector_get_data(borderVector);
600
601 nCols = cpl_image_get_size_x(image);
602 nRows = cpl_image_get_size_y(image);
603 if (nCols < borderSize*2) {
604 BRK_WITH_ERROR_MSG(CPL_ERROR_ILLEGAL_INPUT,
605 "CPL image too small");
606 }
607 type = cpl_image_get_type(image);
608 if (type == CPL_TYPE_FLOAT) {
609 fData = cpl_image_get_data_float(image);
610 } else if (type == CPL_TYPE_DOUBLE) {
611 dData = cpl_image_get_data_double(image);
612 } else {
613 BRK_WITH_ERROR_MSG(CPL_ERROR_INVALID_TYPE,
614 "CPL image type must be float or double");
615 }
616
617 for (cpl_size row=0; row<nRows; row++) {
618 for (int col=0; col<4; col++) {
619 if (type == CPL_TYPE_FLOAT) {
620 borderData[col] = fData[col+row*nCols];
621 borderData[col+borderSize] = fData[nCols-1-col+row*nCols];
622 } else {
623 borderData[col] = dData[col+row*nCols];
624 borderData[col+borderSize] = dData[nCols-1-col+row*nCols];
625 }
626 }
627 borderMean = cpl_vector_get_mean(borderVector);
628 borderMedian = cpl_vector_get_median(borderVector);
629 if (type == CPL_TYPE_FLOAT) {
630 fborderMean = (float) borderMean;
631 fborderMedian = (float) borderMedian;
632 }
633 for (int col=borderSize; col<(nCols-borderSize); col++) {
634 if (type == CPL_TYPE_FLOAT) {
635 fData[col+row*nCols] = fData[col+row*nCols] +
636 fborderMean - fborderMedian;
637 } else {
638 dData[col+row*nCols] = dData[col+row*nCols] +
639 borderMean - borderMedian;
640 }
641 }
642 }
644 }
645 CATCH
646 {
647 retVal = cpl_error_get_code();
648 }
649 eris_ifu_free_vector(&borderVector);
650 eris_check_error_code("eris_ifu_exposure_line_correction");
651 return retVal;
652}
653
654/*----------------------------------------------------------------------------*/
674/*----------------------------------------------------------------------------*/
675cpl_error_code eris_ifu_exposure_column_correction(cpl_image *image)
676{
677 cpl_error_code retVal = CPL_ERROR_NONE;
678 cpl_size nCols;
679 cpl_size nRows;
680 cpl_type type;
681 double *dData = NULL;
682 float *fData= NULL;
683 cpl_vector *borderVector = NULL;
684 double *borderData;
685 double borderMedian;
686 float fborderMedian;
687 const int borderSize = 4;
688
689 cpl_ensure_code(image, CPL_ERROR_NULL_INPUT);
690 TRY
691 {
693 borderVector = cpl_vector_new(borderSize*2));
694 borderData = cpl_vector_get_data(borderVector);
695
696 nCols = cpl_image_get_size_x(image);
697 nRows = cpl_image_get_size_y(image);
698 if (nRows < borderSize*2) {
699 BRK_WITH_ERROR_MSG(CPL_ERROR_ILLEGAL_INPUT,
700 "CPL image too small");
701 }
702 type = cpl_image_get_type(image);
703 if (type == CPL_TYPE_FLOAT) {
704 fData = cpl_image_get_data_float(image);
705 } else if (type == CPL_TYPE_DOUBLE) {
706 dData = cpl_image_get_data_double(image);
707 } else {
708 BRK_WITH_ERROR_MSG(CPL_ERROR_INVALID_TYPE,
709 "CPL image type must be float or double");
710 }
711
712 for (cpl_size col=0; col<nCols; col++) {
713 for (int row=0; row<4; row++) {
714 if (type == CPL_TYPE_FLOAT) {
715 borderData[row] = fData[col + row*nCols];
716 borderData[row+borderSize] = fData[col+(nRows-row-1)*nCols];
717 } else {
718 borderData[row] = dData[col + row*nCols];
719 borderData[row+borderSize] = dData[col+(nRows-row-1)*nCols];
720 }
721 }
722 borderMedian = cpl_vector_get_median(borderVector);
723 if (type == CPL_TYPE_FLOAT) {
724 fborderMedian = (float) borderMedian;
725 }
726 for (int row=borderSize; row<(nCols-borderSize); row++) {
727 if (type == CPL_TYPE_FLOAT) {
728 fData[col + row*nCols] = fData[col + row*nCols] -
729 fborderMedian;
730 } else {
731 dData[col + row*nCols] = dData[col + row*nCols] -
732 borderMedian;
733 }
734 }
735 }
737 }
738 CATCH
739 {
740 retVal = cpl_error_get_code();
741 }
742 eris_ifu_free_vector(&borderVector);
743 eris_check_error_code("eris_ifu_exposure_column_correction");
744 return retVal;
745}
746
747/*----------------------------------------------------------------------------*/
770/*----------------------------------------------------------------------------*/
771cpl_error_code eris_ifu_calc_bpm(const cpl_parameterlist *pl,
772 const char *recipe_name,
773 hdrl_image *master_img,
774 const hdrl_imagelist *imglist_on,
775 cpl_mask **bpm2dMask,
776 cpl_mask **bpm3dMask)
777{
778 cpl_error_code err = CPL_ERROR_NONE;
779 const cpl_parameter *param = NULL;
780 const char *bpmMethod = NULL;
781 char *parName = NULL;
782 hdrl_parameter *bpm2dParam = NULL;
783 hdrl_parameter *bpm3dParam = NULL;
784 cpl_image *bpm3dImg = NULL;
785 cpl_imagelist *bpm3dImageList = NULL;
786 cpl_mask *masterBpmMask = NULL;
787
788 cpl_ensure_code(pl, CPL_ERROR_NULL_INPUT);
789 cpl_ensure_code(recipe_name, CPL_ERROR_NULL_INPUT);
790 cpl_ensure_code(master_img, CPL_ERROR_NULL_INPUT);
791 cpl_ensure_code(imglist_on, CPL_ERROR_NULL_INPUT);
792 cpl_ensure_code(bpm2dMask, CPL_ERROR_NULL_INPUT);
793 cpl_ensure_code(bpm3dMask, CPL_ERROR_NULL_INPUT);
794
795 TRY
796 {
797 parName = cpl_sprintf("eris.%s.bpm.method", recipe_name);
798 param = cpl_parameterlist_find_const(pl, parName);
799 eris_ifu_free_string(&parName);
800
801 bpmMethod = cpl_parameter_get_string(param);
802
803 if (strpbrk(bpmMethod, "2") != NULL)
804 {
805 parName = cpl_sprintf("eris.%s.2dBadPix", recipe_name);
806 bpm2dParam = hdrl_bpm_2d_parameter_parse_parlist(pl, parName);
807 eris_ifu_free_string(&parName);
808
809 cpl_msg_info(cpl_func, "Generating 2D bad pixel mask...");
810 *bpm2dMask = hdrl_bpm_2d_compute(master_img, bpm2dParam);
812 }
813
814 if (strpbrk(bpmMethod, "3") != NULL)
815 {
816 parName = cpl_sprintf("eris.%s.3dBadPix", recipe_name);
817 bpm3dParam = hdrl_bpm_3d_parameter_parse_parlist(pl, parName);
818 eris_ifu_free_string(&parName);
819
820 cpl_msg_info(cpl_func, "Generate 3D bad pixel mask");
821 bpm3dImageList = hdrl_bpm_3d_compute(imglist_on, bpm3dParam);
822
824 int imglistCnt = (int) cpl_imagelist_get_size(bpm3dImageList);
825 cpl_image *patternImg = NULL;
826
827 bpm3dImg = cpl_image_duplicate(cpl_imagelist_get(bpm3dImageList, 0));
828 patternImg = cpl_image_duplicate(
829 cpl_imagelist_get(bpm3dImageList, 0));
830 for (int i = 1; i < imglistCnt; i++) {
831 cpl_image *tmpImg = cpl_imagelist_get(bpm3dImageList, i);
832 cpl_image_add(bpm3dImg, tmpImg);
833 cpl_image_multiply_scalar(tmpImg, (double) (1<<i));
834 cpl_image_add(patternImg, tmpImg);
835 }
836
837 parName = cpl_sprintf("eris.%s.product_depth", recipe_name);
838 param = cpl_parameterlist_find_const(pl, parName);
839 eris_ifu_free_string(&parName);
840 int pd = cpl_parameter_get_int(param);
842 if (pd >= PD_DEBUG) {
843 cpl_image_save(bpm3dImg,ERIS_IFU_PRO_DARK_DBG_3D_FN".fits",
844 CPL_TYPE_INT, NULL, CPL_IO_CREATE);
845 cpl_image_save(patternImg,ERIS_IFU_PRO_DARK_DBG_3D_FN".fits",
846 CPL_TYPE_INT, NULL, CPL_IO_EXTEND);
847 }
848
849 *bpm3dMask = cpl_mask_threshold_image_create(bpm3dImg,
850 1.5, imglistCnt + 0.5);
851 cpl_image_delete(patternImg);
852 }
853
854 cpl_msg_info(cpl_func, "Generating master bad pixel mask...");
855 masterBpmMask = hdrl_image_get_mask(master_img);
856 if (*bpm2dMask != NULL) {
857 cpl_msg_info(cpl_func,"Number of 2D bad pixels: %lld",
858 cpl_mask_count(*bpm2dMask));
859 cpl_mask_or(masterBpmMask, *bpm2dMask);
860 }
861 if (bpm3dImg != NULL) {
862 cpl_msg_info(cpl_func,"Number of 3D bad pixels: %lld",
863 cpl_mask_count(*bpm3dMask));
864 cpl_mask_or(masterBpmMask, *bpm3dMask);
865 }
866 hdrl_image_reject_from_mask(master_img, masterBpmMask);
867 }
868 CATCH
869 {
870 err = cpl_error_get_code();
871
872 cpl_mask_delete(*bpm2dMask); // not using eris_ifu_free_image() here
873 cpl_mask_delete(*bpm3dMask); // because arg is allocated outside of function
874
875 eris_ifu_free_mask(&masterBpmMask);
876 masterBpmMask = cpl_mask_new(hdrl_image_get_size_x(master_img),
877 hdrl_image_get_size_y(master_img));
878 hdrl_image_reject_from_mask(master_img, masterBpmMask);
879 }
880
881 eris_ifu_free_hdrl_parameter(&bpm2dParam);
882 eris_ifu_free_hdrl_parameter(&bpm3dParam);
883 //eris_ifu_free_mask(&masterBpmMask);
884 eris_ifu_free_image(&bpm3dImg);
885 eris_ifu_free_imagelist(&bpm3dImageList);
886 eris_check_error_code("eris_ifu_calc_bpm");
887 return err;
888}
889
890/*----------------------------------------------------------------------------*/
905/*----------------------------------------------------------------------------*/
906hdrl_image *eris_ifu_warp_polynomial_image(const hdrl_image *hdrlInImg,
907 const cpl_polynomial *poly_u,
908 const cpl_polynomial *poly_v)
909{
910 hdrl_image *hdrlOutImg = NULL;
911 const cpl_image *inImgData = NULL,
912 *inImgErr = NULL;
913 cpl_image *outImgData = NULL,
914 *outImgErr = NULL;
915 cpl_vector *profile = NULL;
916 cpl_size nx = 0,
917 ny = 0;
918 cpl_type type = CPL_TYPE_UNSPECIFIED;
919
920 cpl_ensure(hdrlInImg, CPL_ERROR_NULL_INPUT, NULL);
921 cpl_ensure(poly_u, CPL_ERROR_NULL_INPUT, NULL);
922 cpl_ensure(poly_v, CPL_ERROR_NULL_INPUT, NULL);
923
924 TRY
925 {
926 inImgData = hdrl_image_get_image_const(hdrlInImg);
927 inImgErr = hdrl_image_get_error_const(hdrlInImg);
928
929 nx = cpl_image_get_size_x(inImgData);
930 ny = cpl_image_get_size_y(inImgData);
931 type = cpl_image_get_type(inImgData);
933
934 outImgData = cpl_image_new(nx, ny, type);
935 outImgErr = cpl_image_new(nx, ny, type);
936 profile = cpl_vector_new(CPL_KERNEL_DEF_SAMPLES);
937 cpl_vector_fill_kernel_profile(profile, CPL_KERNEL_TANH,
938 CPL_KERNEL_DEF_WIDTH);
939 cpl_image_warp_polynomial(outImgData, inImgData,
940 poly_u, poly_v,
941 profile, CPL_KERNEL_DEF_WIDTH,
942 profile, CPL_KERNEL_DEF_WIDTH);
943 cpl_image_warp_polynomial(outImgErr, inImgErr,
944 poly_u, poly_v,
945 profile, CPL_KERNEL_DEF_WIDTH,
946 profile, CPL_KERNEL_DEF_WIDTH);
947
948 hdrlOutImg = hdrl_image_create(outImgData, outImgErr);
949
950 // there can still be inconsistencies:
951 // warped data can have a value of zero and not be marked as bad
952 // so here we will mark zero data values as bad in the bpm
953 cpl_image *tmpImg = hdrl_image_get_image(hdrlOutImg);
954 double *ptmpImg = cpl_image_get_data(tmpImg);
955 cpl_mask *tmpMask = cpl_image_get_bpm(tmpImg);
956 cpl_binary *ptmpMask = cpl_mask_get_data(tmpMask);
957
958 int w = (int)hdrl_image_get_size_x(hdrlOutImg),
959 h = (int)hdrl_image_get_size_y(hdrlOutImg);
960
961 for (int y = 0; y < h; y++) {
962 for (int x = 0; x < w; x++) {
963 if (ptmpImg[y*w+x] == 0) {
964 ptmpMask[y*w+x] = BAD_PIX;
965 }
966 }
967 }
968 }
969 CATCH
970 {
971 eris_ifu_free_hdrl_image(&hdrlOutImg);
972 }
973 eris_ifu_free_image(&outImgData);
974 eris_ifu_free_image(&outImgErr);
975 eris_ifu_free_vector(&profile);
976 eris_check_error_code("eris_ifu_warp_polynomial_image");
977 return hdrlOutImg;
978}
979
980/*----------------------------------------------------------------------------*/
1001/*----------------------------------------------------------------------------*/
1003 cpl_parameterlist *pl,
1004 const char *recipename)
1005{
1006
1007 cpl_ensure_code(pl,CPL_ERROR_NULL_INPUT);
1008 cpl_ensure_code(recipename,CPL_ERROR_NULL_INPUT);
1009
1010
1011 cpl_error_code retVal = CPL_ERROR_NONE;
1012 cpl_parameter *p = NULL;
1013 char *pName = NULL;
1014 char *recname_full = NULL;
1015 hdrl_parameter *p_hdrl = NULL;
1016 cpl_parameterlist *tmp_pl = NULL;
1017
1018 TRY
1019 {
1020 recname_full = cpl_sprintf("eris.%s", recipename);
1021 pName = cpl_sprintf("%s.%s", recname_full, "instrument");
1022 p = cpl_parameter_new_enum(pName, CPL_TYPE_STRING,
1023 "Specifies the VLT instrument "
1024 "{ERIS,SINFONI,NONE}",
1025 recname_full, "ERIS", 3, "ERIS",
1026 "SINFONI", "NONE");
1027 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "instrument");
1028 cpl_parameterlist_append(pl, p);
1029 p = NULL;
1030 eris_ifu_free_string(&pName);
1031
1032 // productDepthType pdMin = PD_SCIENCE;
1033 // productDepthType pdMax = PD_DEBUG;
1034 pName = cpl_sprintf("%s.%s", recname_full, "product_depth");
1035 p = cpl_parameter_new_value(pName, CPL_TYPE_INT,
1036 "Specifies the product output depth "
1037 "(>0 for auxiliary products)",
1038 recname_full, 0);
1039 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
1040 "product_depth");
1041 cpl_parameterlist_append(pl, p);
1042 p = NULL;
1043 eris_ifu_free_string(&pName);
1044
1045 pName = cpl_sprintf("%s.%s", recname_full, "line_corr");
1046 p = cpl_parameter_new_value(pName, CPL_TYPE_BOOL,
1047 "If TRUE raw exposure image line "
1048 "correction will be applied",
1049 recname_full, CPL_FALSE);
1050 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
1051 "line_corr");
1052 cpl_parameterlist_append(pl, p);
1053 p = NULL;
1054 eris_ifu_free_string(&pName);
1055
1056 pName = cpl_sprintf("%s.%s", recname_full, "col_corr");
1057 p = cpl_parameter_new_value(pName, CPL_TYPE_BOOL,
1058 "If TRUE raw exposure image column "
1059 "correction will be applied",
1060 recname_full, CPL_FALSE);
1061 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
1062 "col_corr");
1063 cpl_parameterlist_append(pl, p);
1064 p = NULL;
1065 eris_ifu_free_string(&pName);
1066
1067 pName = cpl_sprintf("%s.%s", recname_full, "crh_corr");
1068 p = cpl_parameter_new_value(pName, CPL_TYPE_BOOL,
1069 "If TRUE raw exposure image cosmic ray "
1070 "hit correction will be applied",
1071 recname_full, CPL_FALSE);
1072 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
1073 "crh_corr");
1074 cpl_parameterlist_append(pl, p);
1075 p = NULL;
1076 eris_ifu_free_string(&pName);
1077
1078 pName = cpl_sprintf("%s.%s", recname_full, "crh_detection");
1079 p = cpl_parameter_new_value(pName, CPL_TYPE_BOOL,
1080 "If TRUE raw exposure image cosmic ray "
1081 "hit detection will be applied",
1082 recname_full, CPL_FALSE);
1083 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
1084 "crh_detection");
1085 cpl_parameterlist_append(pl, p);
1086 p = NULL;
1087 eris_ifu_free_string(&pName);
1088
1089 p_hdrl = hdrl_lacosmic_parameter_create(5, 2., 3);
1090 tmp_pl = hdrl_lacosmic_parameter_create_parlist(recname_full, "crh",
1091 p_hdrl);
1093
1094
1095 pName = cpl_sprintf("%s.%s", recname_full, "pixel_saturation");
1096 p = cpl_parameter_new_value(pName, CPL_TYPE_DOUBLE,
1097 "Pixel saturation level ",
1098 recname_full, 18000.);
1099 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "pixel_saturation");
1100
1101 cpl_parameterlist_append(pl, p);
1102 p = NULL;
1103
1104
1105
1106 }
1107 CATCH
1108 {
1109 CATCH_MSGS();
1110 retVal = cpl_error_get_code();
1111 }
1112
1116 eris_ifu_free_string(&pName);
1117 eris_ifu_free_string(&recname_full);
1118 eris_check_error_code("eris_ifu_add_std_params");
1119 return retVal;
1120}
1121
1122/*----------------------------------------------------------------------------*/
1131/*----------------------------------------------------------------------------*/
1132void eris_ifu_free_std_param(struct stdParamStruct * stdParams) {
1133
1134 if(stdParams != NULL) {
1135 hdrl_parameter_delete(stdParams->crh_detection);
1136 stdParams=NULL;
1137 }
1138
1139
1140 eris_check_error_code("eris_ifu_free_std_param");
1141 return;
1142}
1143
1144/*----------------------------------------------------------------------------*/
1168/*----------------------------------------------------------------------------*/
1170 const cpl_parameterlist *parlist,
1171 const char *recipename,
1172 struct stdParamStruct *stdParams)
1173{
1174
1175 cpl_ensure_code(parlist,CPL_ERROR_NULL_INPUT);
1176 cpl_ensure_code(recipename,CPL_ERROR_NULL_INPUT);
1177 cpl_ensure_code(stdParams,CPL_ERROR_NULL_INPUT);
1178
1179 cpl_error_code retVal = CPL_ERROR_NONE;
1180 cpl_parameter *p = NULL;
1181 char *pName = NULL;
1182 const char *instrument = NULL;
1183 char *recname_full = NULL;
1184 char *tmp_str = NULL;
1185
1186 TRY
1187 {
1188 recname_full = cpl_sprintf("eris.%s", recipename);
1189 pName = cpl_sprintf("%s.%s", recname_full, "product_depth");
1190 stdParams->productDepth = cpl_parameter_get_int(
1191 cpl_parameterlist_find_const(parlist, pName));
1193 eris_ifu_free_string(&pName);
1194
1195 stdParams->rawImageCorrectionMask = 0;
1196 pName = cpl_sprintf("%s.%s", recname_full, "line_corr");
1197 stdParams->line_corr = cpl_parameter_get_bool(
1198 cpl_parameterlist_find_const(parlist, pName));
1200 eris_ifu_free_string(&pName);
1201 if (stdParams->line_corr) {
1202 stdParams->rawImageCorrectionMask |= LINE_EXPOSURE_CORRECTION;
1203 }
1204
1205 pName = cpl_sprintf("%s.%s", recname_full, "col_corr");
1206 stdParams->col_corr = cpl_parameter_get_bool(
1207 cpl_parameterlist_find_const(parlist, pName));
1209 eris_ifu_free_string(&pName);
1210 if (stdParams->col_corr) {
1211 stdParams->rawImageCorrectionMask |= COLUMN_EXPOSURE_CORRECTION;
1212 }
1213
1214 pName = cpl_sprintf("%s.%s", recname_full, "crh_corr");
1215 stdParams->crh_corr = cpl_parameter_get_bool(
1216 cpl_parameterlist_find_const(parlist, pName));
1218 eris_ifu_free_string(&pName);
1219 if (stdParams->crh_corr) {
1220 stdParams->rawImageCorrectionMask |= COSMIC_RAY_EXPOSURE_CORRECTION;
1221 }
1222
1223 pName = cpl_sprintf("%s.%s", recname_full, "crh_detection");
1224 stdParams->crh_det = cpl_parameter_get_bool(
1225 cpl_parameterlist_find_const(parlist, pName));
1227 eris_ifu_free_string(&pName);
1228 if (stdParams->crh_det) {
1229 stdParams->rawImageCorrectionMask |= COSMIC_RAY_EXPOSURE_DETECTION;
1230 }
1231
1232 tmp_str = cpl_sprintf("%s.crh", recname_full);
1233
1234 stdParams->crh_detection =
1235 hdrl_lacosmic_parameter_parse_parlist(parlist, tmp_str);
1236 eris_ifu_free_string(&tmp_str);
1238
1239 pName = cpl_sprintf("%s.%s", recname_full, "instrument");
1240 instrument = cpl_parameter_get_string(
1241 cpl_parameterlist_find_const(parlist, pName));
1243 eris_ifu_free_string(&pName);
1244 if (strcmp(instrument, "SINFONI") == 0) {
1245 stdParams->instrument = SPIFFI;
1246 }
1247 else if(strcmp(instrument, "ERIS") == 0) {
1248 stdParams->instrument = SPIFFIER;
1249 }
1250 else {
1251 stdParams->instrument = OTHER_INSTRUMENT;
1252 }
1253 }
1254 CATCH
1255 {
1256 CATCH_MSGS();
1257 retVal = cpl_error_get_code();
1258 }
1259
1261 eris_ifu_free_string(&pName);
1262 eris_ifu_free_string(&tmp_str);
1263 eris_ifu_free_string(&recname_full);
1264 eris_check_error_code("eris_ifu_fetch_std_param");
1265 return retVal;
1266}
1267
1268/*----------------------------------------------------------------------------*/
1281/*----------------------------------------------------------------------------*/
1282cpl_error_code eris_parlist_config_add_all_recipes(cpl_parameterlist *pl,
1283 const char* recname)
1284{
1285 char *tmp_str1 = NULL,
1286 *recname_full = NULL;
1287 cpl_error_code err = CPL_ERROR_NONE;
1288 cpl_parameter *p = NULL;
1289
1290 cpl_ensure_code(pl, CPL_ERROR_NULL_INPUT);
1291 cpl_ensure_code(recname, CPL_ERROR_NULL_INPUT);
1292
1293 TRY
1294 {
1295 recname_full = cpl_sprintf("eris.%s", recname);
1296
1297 /* --instrument */
1298 tmp_str1 = cpl_sprintf("%s%s%s", "eris.", recname, ".instrument");
1299 p = cpl_parameter_new_enum(tmp_str1, CPL_TYPE_STRING,
1300 "Specifies the VLT instrument {ERIS,SINFONI,NONE}",
1301 recname_full, INSTRUMENT_ERIS, 3,
1302 INSTRUMENT_ERIS, INSTRUMENT_SINFONI, INSTRUMENT_NONE);
1303 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "instrument");
1304 cpl_parameterlist_append(pl, p);
1305 eris_ifu_free_string(&tmp_str1);
1306
1307 /* -- product_depth */
1308 tmp_str1 = cpl_sprintf("%s%s%s", "eris.", recname, ".product_depth");
1309 p = cpl_parameter_new_value(tmp_str1, CPL_TYPE_INT,
1310 "Specifies the product output depth "
1311 "instrument (>0 for auxillary products)",
1312 recname_full, 0);
1313 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
1314 "product_depth");
1315 cpl_parameterlist_append(pl, p);
1316 eris_ifu_free_string(&tmp_str1);
1317 }
1318 CATCH
1319 {
1320 CATCH_MSGS();
1321 err = cpl_error_get_code();
1323 eris_ifu_free_string(&tmp_str1);
1324 }
1325 eris_ifu_free_string(&tmp_str1);
1326 eris_ifu_free_string(&recname_full);
1327 eris_check_error_code("eris_parlist_config_add_all_recipes");
1328 return err;
1329}
1330
1331/*----------------------------------------------------------------------------*/
1349/*----------------------------------------------------------------------------*/
1350cpl_error_code eris_parlist_config_add_bpm(cpl_parameterlist *pl,
1351 const char* recname)
1352{
1353 char *tmp_str = NULL,
1354 *recname_full = NULL,
1355 *tmp_def_meth = NULL,
1356 *tmp_meth = NULL;
1357 cpl_error_code err = CPL_ERROR_NONE;
1358 cpl_parameterlist *tmp_pl = NULL;
1359 cpl_parameter *p = NULL;
1360 hdrl_parameter *p_hdrl1 = NULL;
1361 hdrl_parameter *p_hdrl2 = NULL;
1362
1363 cpl_ensure_code(pl, CPL_ERROR_NULL_INPUT);
1364 cpl_ensure_code(recname, CPL_ERROR_NULL_INPUT);
1365
1366 TRY
1367 {
1368 recname_full = cpl_sprintf("eris.%s", recname);
1369
1370 /* set default for bpm-method depending on recipe */
1371 if (strcmp(recname, REC_NAME_DISTORTION) == 0) {
1372 tmp_def_meth = cpl_sprintf("2d");
1373 tmp_meth = cpl_sprintf("LEGENDRE");
1374 } else {
1375 tmp_def_meth = cpl_sprintf("2d3d");
1376 tmp_meth = cpl_sprintf("FILTER");
1377 }
1378
1379 /* --bpm.method */
1380 tmp_str = cpl_sprintf("%s%s", recname_full, ".bpm.method");
1381 p = cpl_parameter_new_enum(tmp_str, CPL_TYPE_STRING,
1382 "Specifies the VLT instrument {2d,3d,2d3d}",
1383 recname_full, tmp_def_meth, 4,
1384 "2d", "3d", "2d3d", "none");
1385 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bpm.method");
1386 cpl_parameterlist_append(pl, p);
1387 eris_ifu_free_string(&tmp_str);
1388
1389 /* --collapse.sigclip/minmax (HDRL) */
1390 hdrl_parameter * mode_def =
1391 hdrl_collapse_mode_parameter_create(10., 1., 0., HDRL_MODE_MEDIAN, 0);
1392 p_hdrl1 = hdrl_collapse_sigclip_parameter_create(3., 3., 5);
1393 p_hdrl2 = hdrl_collapse_minmax_parameter_create(1., 1.);
1394 tmp_pl = hdrl_collapse_parameter_create_parlist(recname_full,
1395 "collapse", "MEDIAN",
1396 p_hdrl1, p_hdrl2, mode_def);
1401 hdrl_parameter_delete(mode_def);
1402 /* --2dBadPix.* (HDRL) */
1404 10., 10., 3, CPL_FILTER_MEDIAN, CPL_BORDER_NOP, 5, 5);
1405 // values from HDRL-Demo
1406 // 3, 3, 5, CPL_FILTER_MEDIAN, CPL_BORDER_FILTER, 3, 3));
1408 10., 10., 3, 20, 20, 10, 10, 2, 2);
1409 // values from HDRL-Demo
1410 // 3, 3, 5, 20, 20, 11, 11, 3, 3));
1411 tmp_pl = hdrl_bpm_2d_parameter_create_parlist(recname_full,
1412 "2dBadPix", tmp_meth,
1413 p_hdrl1, p_hdrl2);
1418
1419 /* --3dBadPix.* (HDRL) */
1420 p_hdrl1 = hdrl_bpm_3d_parameter_create(12., 12.,
1421 HDRL_BPM_3D_THRESHOLD_RELATIVE);
1422 tmp_pl = hdrl_bpm_3d_parameter_create_parlist(recname_full,
1423 "3dBadPix", p_hdrl1);
1427 }
1428 CATCH
1429 {
1430 CATCH_MSGS();
1431 err = cpl_error_get_code();
1433 eris_ifu_free_string(&tmp_str);
1437 }
1438 eris_ifu_free_string(&tmp_str);
1439 eris_ifu_free_string(&recname_full);
1440 eris_ifu_free_string(&tmp_def_meth);
1441 eris_ifu_free_string(&tmp_meth);
1442 eris_check_error_code("eris_parlist_config_add_bpm");
1443 return err;
1444}
1445
1446/*----------------------------------------------------------------------------*/
1464/*----------------------------------------------------------------------------*/
1465cpl_error_code eris_parlist_config_add_flat(cpl_parameterlist *pl,
1466 const char* recname)
1467{
1468 cpl_error_code err = CPL_ERROR_NONE;
1469 hdrl_parameter *p_hdrl = NULL;
1470 cpl_parameterlist *tmp_pl = NULL;
1471 cpl_parameter *p = NULL;
1472 char *tmp_str = NULL,
1473 *recname_full = NULL;
1474
1475 cpl_ensure_code(pl, CPL_ERROR_NULL_INPUT);
1476
1477 TRY
1478 {
1479 recname_full = cpl_sprintf("eris.%s", recname);
1480
1481 /* add mode parameter */
1482 tmp_str = cpl_sprintf("%s%s", recname_full, ".mode");
1483 // set default to segment, when this mode is implemented
1484 p = cpl_parameter_new_enum(tmp_str,
1485 CPL_TYPE_STRING,
1486 "Mode of flat-calculation",
1487 recname,
1488 flatModes[0], 3,
1489 flatModes[0],flatModes[1],flatModes[2]);
1490 // if (strcmp(recname, REC_NAME_DISTORTION)==0) {
1491 // BRK_IF_NULL(
1492 // /* fast-mode for distortion */
1493 // p = cpl_parameter_new_enum(tmp_str,
1494 // CPL_TYPE_STRING,
1495 // "Mode of flat-calculation",
1496 // recname,
1497 // flatModes[2], 3,
1498 // flatModes[0],flatModes[1],flatModes[2]));
1499 // } else {
1500 // BRK_IF_NULL(
1501 // /* hdrl-mode for flat */
1502 // p = cpl_parameter_new_enum(tmp_str,
1503 // CPL_TYPE_STRING,
1504 // "Mode of flat-calculation",
1505 // recname,
1506 // flatModes[1], 3,
1507 // flatModes[0],flatModes[1],flatModes[2]));
1508 // }
1509
1510 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,"flat.mode");
1511 cpl_parameterlist_append(pl, p);
1512 eris_ifu_free_string(&tmp_str);
1513
1514 /* --flat_lo.* (HDRL) */
1515 p_hdrl = hdrl_flat_parameter_create(5, 5, HDRL_FLAT_FREQ_LOW);
1516 tmp_pl = hdrl_flat_parameter_create_parlist(recname_full, "flat_lo",
1517 p_hdrl);
1521
1522 /* --flat_hi.* (HDRL) */
1523 p_hdrl = hdrl_flat_parameter_create(7, 7, HDRL_FLAT_FREQ_HIGH);
1524 tmp_pl = hdrl_flat_parameter_create_parlist(recname_full, "flat_hi",
1525 p_hdrl);
1529
1530 /* QC FPN */
1531 tmp_str = cpl_sprintf("%s%s", recname_full, ".qc.fpn.xmin1");
1532 p = cpl_parameter_new_range(tmp_str, CPL_TYPE_INT,
1533 "Fixed Pattern Noise: qc_fpn_xmin1",
1534 recname_full,
1535 512, 1, ERIS_IFU_DETECTOR_SIZE_X);
1536 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "qc_fpn_xmin1");
1537 cpl_parameterlist_append(pl, p);
1538 eris_ifu_free_string(&tmp_str);
1539
1540 tmp_str = cpl_sprintf("%s%s", recname_full, ".qc.fpn.xmax1");
1541 p = cpl_parameter_new_range(tmp_str, CPL_TYPE_INT,
1542 "Fixed Pattern Noise: qc_fpn_xmax1",
1543 recname_full,
1544 1536, 1, ERIS_IFU_DETECTOR_SIZE_X);
1545 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "qc_fpn_xmax1");
1546 cpl_parameterlist_append(pl, p);
1547 eris_ifu_free_string(&tmp_str);
1548
1549 tmp_str = cpl_sprintf("%s%s", recname_full, ".qc.fpn.ymin1");
1550 p = cpl_parameter_new_range(tmp_str, CPL_TYPE_INT,
1551 "Fixed Pattern Noise: qc_fpn_ymin1",
1552 recname_full,
1553 512, 1, ERIS_IFU_DETECTOR_SIZE_Y);
1554 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "qc_fpn_ymin1");
1555 cpl_parameterlist_append(pl, p);
1556 eris_ifu_free_string(&tmp_str);
1557
1558 tmp_str = cpl_sprintf("%s%s", recname_full, ".qc.fpn.ymax1");
1559 p = cpl_parameter_new_range(tmp_str, CPL_TYPE_INT,
1560 "Fixed Pattern Noise: qc_fpn_ymax1",
1561 recname_full,
1562 1536, 1, ERIS_IFU_DETECTOR_SIZE_Y);
1563 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "qc_fpn_ymax1");
1564 cpl_parameterlist_append(pl, p);
1565 eris_ifu_free_string(&tmp_str);
1566
1567 tmp_str = cpl_sprintf("%s%s", recname_full, ".qc.fpn.xmin2");
1568 p = cpl_parameter_new_range(tmp_str, CPL_TYPE_INT,
1569 "Fixed Pattern Noise: qc_fpn_xmin2",
1570 recname_full,
1571 1350, 1, ERIS_IFU_DETECTOR_SIZE_X);
1572 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "qc_fpn_xmin2");
1573 cpl_parameterlist_append(pl, p);
1574 eris_ifu_free_string(&tmp_str);
1575
1576 tmp_str = cpl_sprintf("%s%s", recname_full, ".qc.fpn.xmax2");
1577 p = cpl_parameter_new_range(tmp_str, CPL_TYPE_INT,
1578 "Fixed Pattern Noise: qc_fpn_xmax2",
1579 recname_full,
1580 1390, 1, ERIS_IFU_DETECTOR_SIZE_X);
1581 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "qc_fpn_xmax2");
1582 cpl_parameterlist_append(pl, p);
1583 eris_ifu_free_string(&tmp_str);
1584
1585 tmp_str = cpl_sprintf("%s%s", recname_full, ".qc.fpn.ymin2");
1586 p = cpl_parameter_new_range(tmp_str, CPL_TYPE_INT,
1587 "Fixed Pattern Noise: qc_fpn_ymin2",
1588 recname_full,
1589 1000, 1, ERIS_IFU_DETECTOR_SIZE_Y);
1590 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "qc_fpn_ymin2");
1591 cpl_parameterlist_append(pl, p);
1592 eris_ifu_free_string(&tmp_str);
1593
1594 tmp_str = cpl_sprintf("%s%s", recname_full, ".qc.fpn.ymax2");
1595 p = cpl_parameter_new_range(tmp_str, CPL_TYPE_INT,
1596 "Fixed Pattern Noise: qc_fpn_ymax2",
1597 recname_full,
1598 1200, 1, ERIS_IFU_DETECTOR_SIZE_Y);
1599 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "qc_fpn_ymax2");
1600 cpl_parameterlist_append(pl, p);
1601 eris_ifu_free_string(&tmp_str);
1602 }
1603 CATCH
1604 {
1605 CATCH_MSGS();
1606 err = cpl_error_get_code();
1609 }
1610 eris_ifu_free_string(&tmp_str);
1611 eris_ifu_free_string(&recname_full);
1612 eris_check_error_code("eris_parlist_config_add_flat");
1613 return err;
1614}
1615
1616/*----------------------------------------------------------------------------*/
1631/*----------------------------------------------------------------------------*/
1632cpl_vector* eris_ifu_polyfit_1d(const cpl_vector *x,
1633 const cpl_vector *y,
1634 const int degree)
1635{
1636 cpl_vector *fit_par = NULL;
1637 cpl_polynomial *poly_coeff = NULL;
1638 cpl_matrix *x_matrix = NULL;
1639 double *pfit_par = NULL,
1640 *px = NULL;
1641 cpl_size k = 0,
1642 mindeg1d = 0, //1,
1643 maxdeg1d = degree;
1644
1645 const cpl_boolean sampsym = CPL_FALSE;
1646
1647 cpl_ensure(x, CPL_ERROR_NULL_INPUT, NULL);
1648 cpl_ensure(y, CPL_ERROR_NULL_INPUT, NULL);
1649 cpl_ensure(degree > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
1650
1651 TRY
1652 {
1653
1654 //
1655 // setup data for fitting
1656 //
1657 poly_coeff = cpl_polynomial_new(1);
1658
1659 // put x-vector into a matrix (cast away constness of x: is okay since
1660 // data is wrapped into a matrix which is passed as const again)
1661 px = cpl_vector_get_data((cpl_vector*)x);
1662 x_matrix = cpl_matrix_wrap(1, cpl_vector_get_size(x), px);
1663
1664 //
1665 // fit 1d data
1666 //
1667 cpl_polynomial_fit(poly_coeff,
1668 x_matrix,
1669 &sampsym,
1670 y,
1671 NULL,
1672 CPL_FALSE,
1673 &mindeg1d,
1674 &maxdeg1d);
1675
1676 cpl_matrix_unwrap(x_matrix); x_matrix = NULL;
1678
1679 //
1680 // put fit coefficients into a vector to return
1681 //
1682 fit_par = cpl_vector_new(degree + 1);
1683 pfit_par = cpl_vector_get_data(fit_par);
1684
1685 for(k = 0; k <= degree; k++) {
1686 pfit_par[k] = cpl_polynomial_get_coeff(poly_coeff, &k);
1687 }
1688 }
1689 CATCH
1690 {
1691 eris_ifu_free_vector(&fit_par);
1692 }
1693
1694 cpl_matrix_unwrap(x_matrix); x_matrix = NULL;
1695 cpl_polynomial_delete(poly_coeff); poly_coeff = NULL;
1696 eris_check_error_code("eris_ifu_polyfit_1d");
1697 return fit_par;
1698}
1699
1700/*----------------------------------------------------------------------------*/
1715/*----------------------------------------------------------------------------*/
1717 int nPoints,
1718 double *xdata,
1719 double *ydata,
1720 int degree)
1721{
1722 cpl_polynomial *fit = NULL;
1723 cpl_matrix *x_pos = NULL;
1724 cpl_vector *values = NULL;
1725 TRY
1726 {
1727 fit = cpl_polynomial_new(1);
1728 if (nPoints < degree+1) {
1729 for (cpl_size power=0; power<=degree; power++) {
1730 cpl_polynomial_set_coeff(fit, &power, 0.0);
1731 }
1732
1733 } else {
1734 x_pos = cpl_matrix_wrap(1, nPoints, xdata);
1735 values = cpl_vector_wrap(nPoints, ydata);
1736 const cpl_boolean sampsym = CPL_FALSE;
1737 const cpl_size maxdeg1d = degree;
1738 cpl_polynomial_fit(fit, x_pos, &sampsym, values, NULL,
1739 CPL_FALSE, NULL, &maxdeg1d);
1740 }
1741 }
1742 CATCH
1743 {
1744 CATCH_MSGS();
1745 if (fit != NULL) {
1746 cpl_polynomial_delete(fit);
1747 }
1748 fit = NULL;
1749 }
1750
1751 if (x_pos != NULL) {
1752 cpl_matrix_unwrap(x_pos);
1753 }
1754 if (values != NULL) {
1755 cpl_vector_unwrap(values);
1756 }
1757 eris_check_error_code("eris_ifu_1d_polynomial_fit");
1758 return fit;
1759}
1760
1761/*----------------------------------------------------------------------------*/
1787/*----------------------------------------------------------------------------*/
1789 double *xIn,
1790 double *yIn,
1791 int nIn,
1792 double *xOut,
1793 double *yOut,
1794 int nOut,
1795 const int interType)
1796{
1797 cpl_error_code retVal = CPL_ERROR_NONE;
1798 const gsl_interp_type *interp;
1799 int gsl_status;
1800 gsl_error_handler_t * gslErrorHandler;
1801 double *xIn2 = NULL;
1802 double *yIn2 = NULL;
1803 int nIn2 = 0;
1804
1805 TRY
1806 {
1807 gslErrorHandler = gsl_set_error_handler_off();
1808
1809 ASSURE((xIn != NULL && yIn != NULL && xOut != NULL &&
1810 yOut != NULL),
1811 CPL_ERROR_NULL_INPUT, "no NULL pointer allowed for arrays");
1812 ASSURE((nIn > 0 && nOut > 0),
1813 CPL_ERROR_ILLEGAL_INPUT, "invalid array size input");
1814
1815 remove_2nans(nIn, xIn, yIn, &nIn2, &xIn2, &yIn2);
1816
1817 gsl_interp_accel *acc = gsl_interp_accel_alloc();
1818 if (interType > 2){
1819 const int nPoints = interType;
1820 interp = gsl_interp_polynomial;
1821 gsl_interp *workspace = gsl_interp_alloc(interp, nPoints);
1822 size_t sx;
1823 for (int rx = 0; rx < nOut; rx++) {
1824 sx = gsl_interp_accel_find(acc, xIn2, nIn2,
1825 xOut[rx]);
1826 if (sx > (size_t) (nPoints / 2)) {
1827 sx -= nPoints / 2;
1828 }
1829 if (((cpl_size) sx + nPoints) > nIn2) {
1830 sx = nIn2 - nPoints;
1831 }
1832 gsl_status = gsl_interp_init(workspace, &xIn2[sx], &yIn2[sx],
1833 nPoints);
1834 if (gsl_status != 0) {
1835 BRK_WITH_ERROR_MSG(CPL_ERROR_UNSPECIFIED,
1836 "GSL interpolation routine returned error %d: %s",
1837 gsl_status, gsl_strerror(gsl_status));
1838 }
1839 gsl_status = gsl_interp_eval_e(workspace, &xIn2[sx],
1840 &yIn2[sx], xOut[rx], acc, &yOut[rx]);
1841 if (gsl_status != 0 && gsl_status != GSL_EDOM) {
1842 BRK_WITH_ERROR_MSG(CPL_ERROR_UNSPECIFIED,
1843 "GSL interpolation routine returned error %d: %s",
1844 gsl_status, gsl_strerror(gsl_status));
1845 }
1846 }
1847 gsl_interp_free(workspace);
1848 } else {
1849 switch (interType) {
1850 case 2:
1851 interp = gsl_interp_linear;
1852 break;
1853 case -1:
1854 interp = gsl_interp_cspline;
1855 break;
1856 case -2:
1857 interp = gsl_interp_cspline_periodic;
1858 break;
1859 case -3:
1860 interp = gsl_interp_akima;
1861 break;
1862 case -4:
1863 interp = gsl_interp_akima_periodic;
1864 break;
1865 case -5:
1866 interp = gsl_interp_steffen;
1867 break;
1868 default:
1869 BRK_WITH_ERROR_MSG(CPL_ERROR_UNSUPPORTED_MODE,
1870 "unknown interpolation type %d", interType);
1871 }
1872 gsl_interp *workspace = gsl_interp_alloc(interp, nIn2);
1873
1874 gsl_status = gsl_interp_init(workspace, xIn2, yIn2, nIn2);
1875 if (gsl_status != 0) {
1876 BRK_WITH_ERROR_MSG(CPL_ERROR_UNSPECIFIED,
1877 "GSL interpolation routine returned error %d: %s",
1878 gsl_status, gsl_strerror(gsl_status));
1879 }
1880 for (int rx = 0; rx < nOut; rx++) {
1881 gsl_status = gsl_interp_eval_e(workspace, xIn2, yIn2,
1882 xOut[rx], acc, &yOut[rx]);
1883 if (gsl_status != 0 && gsl_status != GSL_EDOM) {
1884 BRK_WITH_ERROR_MSG(CPL_ERROR_UNSPECIFIED,
1885 "GSL interpolation routine returned error %d: %s",
1886 gsl_status, gsl_strerror(gsl_status));
1887 }
1888 }
1889 gsl_interp_free(workspace);
1890 }
1891 gsl_interp_accel_free(acc);
1892
1893 }
1894 CATCH
1895 {
1896 retVal = cpl_error_get_code();
1897 }
1898 gsl_set_error_handler(gslErrorHandler);
1901 eris_check_error_code("eris_ifu_1d_interpolation");
1902 return retVal;
1903}
1904
1905void remove_2nans(int size_in, double *xin, double *yin, int *size_out, double **xout, double **yout) {
1906 int no_valids = 0,
1907 ox = 0;
1908
1909 for (int i = 0; i < size_in; i++) {
1910 if ((! eris_ifu_is_nan_or_inf(xin[i])) && (! eris_ifu_is_nan_or_inf(yin[i]))) {
1911 // if ((kmclipm_is_nan_or_inf(xin[i]) == 0) && (kmclipm_is_nan_or_inf(yin[i]) == 0)) {
1912 no_valids++;
1913 }
1914 }
1915
1916 *size_out = no_valids;
1917 *xout = (double *) cpl_calloc(no_valids, sizeof(double));
1918 *yout = (double *) cpl_calloc(no_valids, sizeof(double));
1919
1920 ox = 0;
1921 for (int i = 0; i < size_in; i++) {
1922 if ((! eris_ifu_is_nan_or_inf(xin[i])) && (! eris_ifu_is_nan_or_inf(yin[i]))) {
1923 // if ((kmclipm_is_nan_or_inf(xin[i]) == 0) && (kmclipm_is_nan_or_inf(yin[i]) == 0)) {
1924 (*xout)[ox] = xin[i];
1925 (*yout)[ox] = yin[i];
1926 ox++;
1927 }
1928 }
1929 eris_check_error_code("remove_2nans");
1930 return;
1931}
1932
1933
1934// now in eris_ifu_vector.h/.c
1935//
1937// @brief Compute the mean value of vector elements ignoring NaN values.
1938// @param v Input const cpl_vector
1939// @return Mean value of vector elements or undefined on error.
1940// */
1941//double eris_ifu_vector_get_mean(const cpl_vector *v)
1942//{
1943// double sum = 0.,
1944// mean = 0.;
1945// const double *pv = NULL;
1946// int n = 0;
1947
1948// cpl_ensure(v != NULL, CPL_ERROR_NULL_INPUT, -1.0);
1949
1950// TRY
1951// {
1952// BRK_IF_NULL(
1953// pv = cpl_vector_get_data_const(v));
1954
1955// for (int i = 0; i < cpl_vector_get_size(v); i++) {
1956// if (!isnan(pv[i])) {
1957// sum += pv[i];
1958// n++;
1959// }
1960// }
1961// mean = sum / (double)n;
1962// }
1963// CATCH
1964// {
1965// CATCH_MSGS();
1966// mean = 0./0.;
1967// }
1968// return mean;
1969//}
1970
1979double eris_ifu_image_get_mean(const cpl_image *image)
1980{
1981 cpl_size nr_valid_pix = 0,
1982 nx = 0,
1983 ny = 0;
1984 const double *pimg = NULL;
1985 double mean = 0.0,
1986 sum = 0.0;
1987
1988 cpl_ensure(image != NULL, CPL_ERROR_NULL_INPUT, -1.0);
1989
1990 TRY
1991 {
1992 nx = cpl_image_get_size_x(image);
1993 ny = cpl_image_get_size_y(image);
1994 nr_valid_pix = nx * ny;
1995
1997 pimg = cpl_image_get_data_double_const(image));
1998
1999 for (int j = 0; j < ny; j++) {
2000 for (int i = 0; i < nx; i++) {
2001 if (isnan(pimg[i+j*nx]))
2002 nr_valid_pix--;
2003 else
2004 sum += pimg[i+j*nx];
2005 }
2006 }
2007 mean = sum / (double) nr_valid_pix;
2008 }
2009 CATCH
2010 {
2011 CATCH_MSGS();
2012 mean = 0.0;
2013 }
2014 eris_check_error_code("eris_ifu_image_get_mean");
2015 return mean;
2016}
2017
2018cpl_error_code eris_ifu_line_gauss_fit(
2019 const cpl_vector *yIn,
2020 int center,
2021 int range,
2022 struct gaussParStruct *gaussPar)
2023{
2024 cpl_error_code retVal = CPL_ERROR_NONE;
2025 int start = 0;
2026 double pos = 0.,
2027 sigma = 0.,
2028 area = 0.,
2029 offset = 0.,
2030 mserr = 0.,
2031 *yInData = NULL,
2032 *pxVec = NULL;
2033 const double *pyVec = NULL;
2034 cpl_vector *xVec = NULL,
2035 *yVec = NULL,
2036 *yIn2 = NULL;
2037
2038 cpl_ensure_code(yIn, CPL_ERROR_NULL_INPUT);
2039 cpl_ensure_code(cpl_vector_get_size(yIn) == ERIS_IFU_DETECTOR_SIZE,
2040 CPL_ERROR_ILLEGAL_INPUT);
2041
2042 TRY
2043 {
2044 try_again_fit:
2046 yIn2 = cpl_vector_duplicate(yIn));
2048 yInData = cpl_vector_get_data(yIn2));
2049
2050 // calc start position in vector-data
2051 start = center-(range/2);
2052
2053 if (start < 0) {
2054 // no negative starting points!
2055 start = 0;
2056 range -= start;
2057 } else if (start >= ERIS_IFU_DETECTOR_SIZE_X-ERIS_IFU_DETECTOR_BP_BORDER) {
2058 start = ERIS_IFU_DETECTOR_SIZE_X - ERIS_IFU_DETECTOR_BP_BORDER - 1;
2059 }
2060 if (start + range >= ERIS_IFU_DETECTOR_SIZE_X-ERIS_IFU_DETECTOR_BP_BORDER) {
2061 // start position + range to fit are larger than input data!
2062 // decrease range accordingly
2063 range = ERIS_IFU_DETECTOR_SIZE_X-ERIS_IFU_DETECTOR_BP_BORDER - start;
2064 }
2066 yVec = cpl_vector_wrap(range, &yInData[start]));
2067
2069 pyVec = cpl_vector_get_data_const(yVec));
2070
2072 xVec = cpl_vector_new(range));
2074 pxVec = cpl_vector_get_data(xVec));
2075 for (int ix=0; ix<range; ix++) {
2076 pxVec[ix] = start + ix;
2077 }
2078
2079 cpl_errorstate prestate = cpl_errorstate_get();
2080 retVal = cpl_vector_fit_gaussian(xVec, NULL, yVec, NULL, CPL_FIT_ALL,
2081 &pos, &sigma, &area, &offset, &mserr, NULL, NULL);
2082
2083 gaussPar->errorCode = (int) retVal;
2084 if ((retVal == CPL_ERROR_NONE) ||
2085 ((retVal == CPL_ERROR_CONTINUE) &&
2086 (((pos > 2042) && (pos < ERIS_IFU_DETECTOR_SIZE_X-ERIS_IFU_DETECTOR_BP_BORDER)) ||
2087 ((pos < 7) && (pos >= ERIS_IFU_DETECTOR_BP_BORDER)))))
2088 {
2089 // retVal is perhaps bad, but fit is close enough for continuing
2090 if (retVal == CPL_ERROR_CONTINUE) {
2091 cpl_errorstate_set(prestate);
2092 retVal = CPL_ERROR_NONE;
2093 }
2094
2095 gaussPar->x0 = pos;
2096 gaussPar->sigma = sigma;
2097 gaussPar->area = area;
2098 gaussPar->offset = offset;
2099 gaussPar->peak = area / sqrt(2 * CPL_MATH_PI * sigma * sigma);
2100 gaussPar->mse = mserr;
2101
2102 } else {
2103 gaussPar->x0 = 0.;
2104 gaussPar->sigma = 0.;
2105 gaussPar->area = 0.;
2106 gaussPar->offset = 0.;
2107 gaussPar->peak = 0.;
2108 gaussPar->mse = 0.;
2109 }
2110
2111 if (retVal == CPL_ERROR_CONTINUE) {
2112 // Recover from "The iterative process did not converge" error
2113 cpl_errorstate_set(prestate);
2114 retVal = CPL_ERROR_NONE;
2115
2116 const int decrement = 9;
2117 if ((range - decrement) > 2) {
2118 range -= decrement; // decrease range,
2119 // in doubt that some artifact is disturbing the fit
2120 cpl_vector_unwrap(yVec); yVec = NULL;
2121 eris_ifu_free_vector(&xVec);
2122 eris_ifu_free_vector(&yIn2);
2123
2124 goto try_again_fit;
2125 }
2126 }
2127
2128 if (range <= GAUSS_PAR_RANGE_MAX) {
2129 for (int ix=0; ix<range; ix++) {
2130 gaussPar->range = range;
2131 gaussPar->xdata[ix] = pxVec[ix];
2132 gaussPar->ydata[ix] = pyVec[ix];
2133 }
2134 }
2135 }
2136 CATCH
2137 {
2138 gaussPar->x0 = 0.;
2139 gaussPar->sigma = 0.;
2140 gaussPar->area = 0.;
2141 gaussPar->offset = 0.;
2142 gaussPar->peak = 0.;
2143 gaussPar->mse = 0.;
2144 if (range <= GAUSS_PAR_RANGE_MAX) {
2145 for (int ix=0; ix<range; ix++) {
2146 gaussPar->range = 0.;
2147 gaussPar->xdata[ix] = 0.;
2148 gaussPar->ydata[ix] = 0.;
2149 }
2150 }
2151 }
2152 eris_ifu_free_vector(&xVec);
2153 eris_ifu_free_vector(&yIn2);
2154 if (yVec != NULL) {
2155 cpl_vector_unwrap(yVec); yVec = NULL;
2156 }
2157 eris_check_error_code("eris_ifu_line_gauss_fit");
2158 return retVal;
2159}
2160
2164cpl_vector* eris_ifu_calc_centers_collapse_chunk(const cpl_image* img,
2165 int chunk_center,
2166 int height)
2167{
2168 cpl_size nx = 0,
2169 ny = 0,
2170 chunk_top = 0,
2171 chunk_bottom = 0;
2172 double *pcollapsed_img = NULL;
2173 cpl_vector *tmp = NULL,
2174 *collapsed = NULL;
2175 cpl_image *collapsed_img = NULL;
2176
2177 cpl_ensure(chunk_center > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
2178 cpl_ensure(height > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
2179 // height should be even
2180 cpl_ensure(height % 2 == 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
2181 cpl_ensure(img, CPL_ERROR_NULL_INPUT, NULL);
2182
2183 TRY
2184 {
2185 nx = cpl_image_get_size_x(img);
2186 ny = cpl_image_get_size_y(img);
2187
2188 /* lower and upper limit of chunk in y */
2189 chunk_bottom = chunk_center - height/2;
2190 chunk_top = ny - (chunk_center + height/2);
2191
2192 /* collapse image in y and convert to vector */
2193 collapsed_img = cpl_image_collapse_median_create(img, 0,
2194 chunk_bottom,
2195 chunk_top);
2196 pcollapsed_img = cpl_image_get_data_double(collapsed_img);
2197 tmp = cpl_vector_wrap(nx, pcollapsed_img);
2198 collapsed = cpl_vector_duplicate(tmp);
2199 }
2200 CATCH
2201 {
2202 CATCH_MSGS();
2203 eris_ifu_free_vector(&collapsed);
2204 }
2205
2206 cpl_vector_unwrap(tmp); tmp = NULL;
2207 eris_ifu_free_image(&collapsed_img);
2208 eris_check_error_code("eris_ifu_calc_centers_collapse_chunk");
2209 return collapsed;
2210}
2211
2215cpl_vector* eris_ifu_image_collapse(const cpl_image* img)
2216{
2217 cpl_size nx = 0;
2218 double *pcollapsed_img = NULL;
2219 cpl_vector *tmp = NULL,
2220 *collapsed = NULL;
2221 cpl_image *collapsed_img = NULL;
2222
2223 cpl_ensure(img, CPL_ERROR_NULL_INPUT, NULL);
2224
2225 TRY
2226 {
2227 nx = cpl_image_get_size_x(img);
2228
2229 /* collapse image in y and convert to vector */
2230 collapsed_img = cpl_image_collapse_median_create(img, 0, 0, 0);
2231 pcollapsed_img = cpl_image_get_data_double(collapsed_img);
2232 tmp = cpl_vector_wrap(nx, pcollapsed_img);
2233 collapsed = cpl_vector_duplicate(tmp);
2234 }
2235 CATCH
2236 {
2237 CATCH_MSGS();
2238 eris_ifu_free_vector(&collapsed);
2239 }
2240
2241 cpl_vector_unwrap(tmp); tmp = NULL;
2242 eris_ifu_free_image(&collapsed_img);
2243 eris_check_error_code("eris_ifu_image_collapse");
2244 return collapsed;
2245}
2246
2260cpl_error_code eris_ifu_slitpos_gauss(const cpl_image *profile_x,
2261 double *left_edge_pos,
2262 double *right_edge_pos,
2263 int llx,
2264 int productDepth)
2265{
2266 cpl_error_code err = CPL_ERROR_NONE;
2267 cpl_vector *x_left = NULL,
2268 *x_right = NULL,
2269 *y_left = NULL,
2270 *y_right = NULL,
2271 *y_left2 = NULL,
2272 *y_right2 = NULL,
2273 *profile_x_v = NULL,
2274 *tmpVector = NULL;
2275 int s = 0,
2276 s2 = 0;
2277 double *px_left = NULL,
2278 *px_right = NULL,
2279 *py_left = NULL,
2280 *py_right = NULL,
2281 x0 = 0.,
2282 sigma = 0.,
2283 area = 0.,
2284 offset = 0.,
2285 median = 0.;
2286 const double *pprofile_x = NULL;
2287 cpl_size maxpos = 0,
2288 minpos = 0,
2289 inBetween = 0;
2290
2291 cpl_ensure_code(profile_x, CPL_ERROR_NULL_INPUT);
2292 cpl_ensure_code(left_edge_pos, CPL_ERROR_NULL_INPUT);
2293 cpl_ensure_code(right_edge_pos, CPL_ERROR_NULL_INPUT);
2294
2295 TRY
2296 {
2297 s = (int) cpl_image_get_size_x(profile_x);
2298 s2 = (int)(s/2);
2299
2300 // create x-Vector (profile_x goes into y :-)
2301 // one for the left edge with 0-47
2302 // one for the right edge with 47-95
2303 // split profile_x into two and convert to vector
2304 x_left = cpl_vector_new(s2);
2305 x_right = cpl_vector_new(s2);
2306 y_left = cpl_vector_new(s2);
2307 y_right = cpl_vector_new(s2);
2308
2309 px_left = cpl_vector_get_data(x_left);
2310 px_right = cpl_vector_get_data(x_right);
2311 py_left = cpl_vector_get_data(y_left);
2312 py_right = cpl_vector_get_data(y_right);
2313
2314 //erw
2315 // BRK_IF_NULL(
2316 // pprofile_x = cpl_image_get_data_double_const(profile_x));
2317 tmpVector = cpl_vector_new_from_image_row(profile_x,1);
2318 profile_x_v = cpl_vector_filter_lowpass_create(
2319 tmpVector, CPL_LOWPASS_LINEAR, 2);
2320 eris_ifu_free_vector(&tmpVector);
2321 pprofile_x = cpl_vector_get_data(profile_x_v);
2322
2323 for (int i = 0; i < s2; i++) {
2324 px_left[i] = i;
2325 // subtract left value from right value (we need a positive peak)
2326 py_left[i] = pprofile_x[i+1] - pprofile_x[i];
2327 px_right[i] = i+s2+1;
2328 // subtract right value from left value (we need a positive peak)
2329 py_right[i] = pprofile_x[i+s2-1] - pprofile_x[i+s2];
2330 }
2331
2332 y_left2 = cpl_vector_filter_lowpass_create(y_left,
2333 CPL_LOWPASS_GAUSSIAN, 5);
2334 y_right2 = cpl_vector_filter_lowpass_create(y_right,
2335 CPL_LOWPASS_GAUSSIAN, 5);
2336
2337 // check that left differences don't start with a negative peak
2338 if (cpl_vector_get_min(y_left2) < -100) {
2339 maxpos = cpl_vector_get_maxpos(y_left2);
2340 cpl_vector *tmpVec = cpl_vector_extract(y_left2,0,maxpos,1);
2341 minpos = cpl_vector_get_minpos(tmpVec);
2342 cpl_vector_delete(tmpVec);
2343 median = cpl_vector_get_median_const(y_left2);
2344 if (cpl_vector_get(y_left2,minpos < median-10.)) {
2345 // if (maxpos > minpos) {
2346 // negative peak is to the left --> arc lamp from other slitlet
2347 inBetween = (cpl_size) (((double) (maxpos + minpos) / 2.) + .5);
2348 for (cpl_size vx = 0; vx < inBetween; vx++) {
2349 if (cpl_vector_get(y_left2, vx) < -10.) {
2350 cpl_vector_set(y_left2, vx, median);
2351 }
2352 }
2353 minpos = cpl_vector_get_minpos(y_left2);
2354 if (cpl_vector_get(y_left2,minpos < median-10.)) {
2355 inBetween = (cpl_size) (((double)(maxpos + minpos) / 2.) + .5);
2356 for (cpl_size vx = inBetween;
2357 vx < cpl_vector_get_size(y_left2); vx++) {
2358 if (cpl_vector_get(y_left2, vx) < -10.) {
2359 cpl_vector_set(y_left2, vx, median);
2360 }
2361 }
2362
2363 }
2364
2365 //BRK_WITH_ERROR(CPL_ERROR_EOL);
2366 } else {
2367 cpl_size vSize = cpl_vector_get_size(y_left2);
2368 for (cpl_size vx = maxpos; vx < vSize; vx++) {
2369 if (cpl_vector_get(y_left2, vx) < -10.) {
2370 cpl_vector_set(y_left2, vx, median);
2371 }
2372 }
2373 }
2374 }
2376
2377 // check that right differences don't end with a negative peak
2378 if (cpl_vector_get_min(y_right2) < -100.) {
2379 maxpos = cpl_vector_get_maxpos(y_right2);
2380 minpos = cpl_vector_get_minpos(y_right2);
2381 median = cpl_vector_get_median_const(y_right2);
2382 if (maxpos < minpos) {
2383 // negative peak is to the right--> arc lamp from other slitlet
2384 cpl_size vSize = cpl_vector_get_size(y_right2);
2385 inBetween = (cpl_size) (((double) (maxpos + minpos) / 2.) + .5);
2386 for (cpl_size vx = inBetween; vx < vSize; vx++) {
2387 if (cpl_vector_get(y_left2, vx) < -10.) {
2388 cpl_vector_set(y_left2, vx, median);
2389 }
2390 }
2391 BRK_WITH_ERROR(CPL_ERROR_EOL);
2392 } else {
2393 for (cpl_size vx=maxpos; vx<=0; vx--) {
2394 if (cpl_vector_get(y_right2, vx) < -10.) {
2395 cpl_vector_set(y_right2, vx, median);
2396 }
2397 }
2398 }
2399 }
2401
2402 if (productDepth & 8) {
2403 cpl_propertylist *pl = NULL;
2405 "eris_ifu_distortion_dbg_slitpos_profile_left_xy.fits",
2406 CPL_IO_CREATE, pl);
2408 "eris_ifu_distortion_dbg_slitpos_profile_left_xy.fits",
2409 CPL_IO_EXTEND, pl);
2411 "eris_ifu_distortion_dbg_slitpos_profile_right_xy.fits",
2412 CPL_IO_CREATE, pl);
2414 "eris_ifu_distortion_dbg_slitpos_profile_right_xy.fits",
2415 CPL_IO_EXTEND, pl);
2417 "eris_ifu_distortion_dbg_slitpos_profile_left_xy.fits",
2418 CPL_IO_EXTEND, pl);
2419 eris_ifu_save_vector_dbg(y_right2,
2420 "eris_ifu_distortion_dbg_slitpos_profile_right_xy.fits",
2421 CPL_IO_EXTEND, pl);
2422 }
2423
2424
2425 // fit gauss-function to left edge
2426 // cpl_msg_debug(cpl_func, " === FIT LEFT EDGE ======"
2427 // "=======================");
2428 err = eris_ifu_fit_gauss(x_left, y_left2, &x0, &sigma, &area, &offset);
2429 if (err != CPL_ERROR_NONE) {
2430 eris_ifu_free_vector(&profile_x_v);
2431 eris_ifu_free_vector(&tmpVector);
2432 eris_ifu_free_vector(&x_left);
2433 eris_ifu_free_vector(&x_right);
2434 eris_ifu_free_vector(&y_left);
2435 eris_ifu_free_vector(&y_right);
2436 eris_ifu_free_vector(&y_left2);
2437 eris_ifu_free_vector(&y_right2);
2438 return err;
2439 }
2440
2441 if (sigma > 25) {
2442 BRK_WITH_ERROR(CPL_ERROR_EOL);
2443 }
2444 // correct half pixel because of subtraction of y-values at beginning
2445 *left_edge_pos = llx + x0 + 0.5;
2446
2447 // fit gauss-function to right edge
2448 // cpl_msg_debug(cpl_func, " === FIT RIGHT EDGE ============================");
2449 err = eris_ifu_fit_gauss(x_right, y_right2, &x0, &sigma, &area, &offset);
2450 if (err != CPL_ERROR_NONE) {
2451 eris_ifu_free_vector(&profile_x_v);
2452 eris_ifu_free_vector(&tmpVector);
2453 eris_ifu_free_vector(&x_left);
2454 eris_ifu_free_vector(&x_right);
2455 eris_ifu_free_vector(&y_left);
2456 eris_ifu_free_vector(&y_right);
2457 eris_ifu_free_vector(&y_left2);
2458 eris_ifu_free_vector(&y_right2);
2459 return err;
2460 }
2461
2462 if (sigma > 25) {
2463 BRK_WITH_ERROR(CPL_ERROR_EOL);
2464 }
2465 // correct half pixel because of subtraction of y-values at beginning
2466 *right_edge_pos = llx + x0 - 0.5;
2467
2468 // cpl_msg_debug(cpl_func, "left/right edge : %g / %g",
2469 // *left_edge_pos, *right_edge_pos);
2470 }
2471 CATCH
2472 {
2473 err = cpl_error_get_code();
2474 if (err != CPL_ERROR_EOL) {
2475 CATCH_MSGS();
2476 }
2477 }
2478 eris_ifu_free_vector(&profile_x_v);
2479 eris_ifu_free_vector(&tmpVector);
2480 eris_ifu_free_vector(&x_left);
2481 eris_ifu_free_vector(&x_right);
2482 eris_ifu_free_vector(&y_left);
2483 eris_ifu_free_vector(&y_right);
2484 eris_ifu_free_vector(&y_left2);
2485 eris_ifu_free_vector(&y_right2);
2486 eris_check_error_code("eris_ifu_slitpos_gauss");
2487 return err;
2488}
2489
2499cpl_error_code eris_ifu_bpm_correction(hdrl_image *himg,
2500 hdrl_image *badPixelMaskImg)
2501{
2502 cpl_mask *m = cpl_mask_new(3,3); // 3x3 mask for filtering; use a cross pattern
2503 cpl_mask_set(m, 1,2, BAD_PIX);
2504 cpl_mask_set(m, 2,1, BAD_PIX);
2505 cpl_mask_set(m, 2,3, BAD_PIX);
2506 cpl_mask_set(m, 3,2, BAD_PIX);
2507
2508 cpl_image *filterd_image = cpl_image_duplicate(hdrl_image_get_image_const(himg));
2509 cpl_image_filter_mask(filterd_image, hdrl_image_get_image_const(himg), m, CPL_FILTER_AVERAGE,
2510 CPL_BORDER_FILTER);
2511
2512 cpl_image *filterd_error = cpl_image_duplicate(hdrl_image_get_error_const(himg));
2513 cpl_image *variance = cpl_image_duplicate(hdrl_image_get_error_const(himg));
2514 cpl_image_power(variance, 2.0);
2515 cpl_image_filter_mask(filterd_error, variance, m, CPL_FILTER_AVERAGE,
2516 CPL_BORDER_FILTER);
2517 cpl_image_delete(variance);
2518
2519 cpl_image_divide_scalar(filterd_error, (double) cpl_mask_count(m));
2520 cpl_image_power(filterd_error, 0.5);
2521
2522 if (badPixelMaskImg == NULL){
2523 cpl_image *badPixelMaskImg1 = cpl_image_new_from_mask(hdrl_image_get_mask_const(himg));
2524 cpl_image *badPixelMaskImg2 = cpl_image_cast(badPixelMaskImg1, CPL_TYPE_DOUBLE);
2525 cpl_image_delete(badPixelMaskImg1);
2526
2527 badPixelMaskImg = hdrl_image_new(
2529 hdrl_image_get_size_y(himg));
2530 hdrl_image_insert(badPixelMaskImg,
2531 badPixelMaskImg2, NULL, 1, 1);
2532
2533 cpl_image_delete(badPixelMaskImg2);
2534 }
2535
2536 cpl_image *filtered_mask = cpl_image_duplicate(hdrl_image_get_image_const(badPixelMaskImg));
2537
2538 cpl_mask_set(m, 2, 2, GOOD_PIX);
2539 cpl_image_filter_mask(filtered_mask, hdrl_image_get_image_const(badPixelMaskImg), m, CPL_FILTER_AVERAGE,
2540 CPL_BORDER_FILTER);
2541
2542
2543 // Get the places to be corrected
2544 double corr_thesh = 0.4;
2545 cpl_image_threshold(filtered_mask, corr_thesh, corr_thesh, 1.0, 0.0);
2546 cpl_image_multiply(filtered_mask, hdrl_image_get_image_const(badPixelMaskImg)); // To be corrected set to 1
2547 cpl_mask *corr_mask_t = cpl_mask_threshold_image_create(filtered_mask,-0.1, 0.1);
2548 cpl_mask_not(corr_mask_t);
2549 cpl_msg_info(__func__, "No. bad pixels to be corrected: %d", (int)cpl_mask_count(corr_mask_t));
2550 cpl_image_multiply(filterd_image, filtered_mask);
2551 cpl_image_multiply(filterd_error, filtered_mask);
2552
2553 cpl_image_subtract_scalar(filtered_mask, 1);
2554 cpl_image_multiply_scalar(filtered_mask, -1); // To be kept: set to 1;
2555 cpl_mask *corr_mask = cpl_mask_threshold_image_create(filtered_mask,-0.1, 0.1);
2556 cpl_mask_not(corr_mask);
2557 cpl_msg_debug(__func__, "No. bad pixels before correction: %d", (int)cpl_image_count_rejected(hdrl_image_get_image(himg)));
2558 cpl_mask_and(hdrl_image_get_mask(himg), corr_mask);
2559 cpl_mask_and(cpl_image_get_bpm(hdrl_image_get_error(himg)), corr_mask);
2560 cpl_msg_debug(__func__, "No. bad pixels after correction: %d", (int)cpl_image_count_rejected(hdrl_image_get_image(himg)));
2561
2562 cpl_image_multiply(hdrl_image_get_image(himg), filtered_mask);
2563 cpl_image_add(hdrl_image_get_image(himg), filterd_image);
2564
2565 cpl_image_multiply(hdrl_image_get_error(himg), filtered_mask);
2566 cpl_image_add(hdrl_image_get_error(himg), filterd_error);
2567
2568 cpl_image_multiply(hdrl_image_get_image(badPixelMaskImg), filtered_mask);
2569 cpl_mask_delete(m);
2570 cpl_image_delete(filtered_mask);
2571 cpl_image_delete(filterd_error);
2572 cpl_image_delete(filterd_image);
2573 cpl_mask_delete(corr_mask);
2574 cpl_mask_delete(corr_mask_t);
2575 eris_check_error_code("eris_ifu_bpm_correction");
2576 return cpl_error_get_code();
2577}
2578
2579double eris_ifu_auto_derot_corr(double alt, double rot)
2580{
2581 double za = 90.0-alt;
2582 double offset = -0.000379*za + 0.000500*pow(za,2);
2583 double alt_scale = 0.0467*za - 0.000265*pow(za,2);
2584 double pi = 3.14159265359;
2585 double derot_corr = offset + alt_scale * cos((rot - 15.4)*pi/180);
2586 eris_check_error_code("eris_ifu_auto_derot_corr");
2587 return derot_corr;
2588}
ifsInstrument eris_ifu_get_instrument(const cpl_propertylist *header)
Return the used instrument of the FITS file header.
cpl_error_code eris_ifu_fit_gauss(const cpl_vector *x, const cpl_vector *y, double *x0, double *sigma, double *area, double *offset)
Fit Gaussian to find peak center and width.
#define BRK_WITH_ERROR(code)
Set a new CPL error, and exit the try-block.
#define ASSURE(condition, error,...)
error handling macro (from fors-pipeline)
#define BRK_IF_ERROR(function)
If function is or returns an error != CPL_ERROR_NONE, then the try-block is exited.
#define CHECK_ERROR_STATE(void)
Check the CPL error state, and exit the try-block if not CPL_ERROR_NONE.
#define BRK_WITH_ERROR_MSG(code,...)
Set a new CPL error, and exit the try-block.
#define CATCH_MSG()
Displays an error message.
#define TRY
Beginning of a TRY-block.
#define CATCH
End of a TRY-block, beginning of a CATCH-block.
#define BRK_IF_NULL(function)
If function is or returns a NULL pointer, then the try-block is exited.
#define CATCH_MSGS()
Displays an error message stack.
cpl_error_code eris_ifu_exposure_line_correction(cpl_image *image)
Perform line/row correction on raw detector image.
cpl_polynomial * eris_ifu_1d_polynomial_fit(int nPoints, double *xdata, double *ydata, int degree)
Fit a 1D polynomial to arrays of data points.
cpl_error_code eris_ifu_slitpos_gauss(const cpl_image *profile_x, double *left_edge_pos, double *right_edge_pos, int llx, int productDepth)
eris_ifu_dist_slitpos_gauss
cpl_error_code eris_ifu_1d_interpolation(double *xIn, double *yIn, int nIn, double *xOut, double *yOut, int nOut, const int interType)
Perform 1D interpolation using GSL routines.
cpl_error_code eris_ifu_bpm_correction(hdrl_image *himg, hdrl_image *badPixelMaskImg)
eris_ifu_bpm_correction
cpl_error_code eris_parlist_config_add_flat(cpl_parameterlist *pl, const char *recname)
Add flat field configuration parameters to parameter list.
cpl_error_code eris_parlist_config_add_all_recipes(cpl_parameterlist *pl, const char *recname)
Add common configuration parameters for all recipes.
cpl_vector * eris_ifu_calc_centers_collapse_chunk(const cpl_image *img, int chunk_center, int height)
hdrl_image * eris_ifu_load_exposure_frame(const cpl_frame *frame, int exposureCorrectionMode, cpl_image *dqi)
Load a raw detector exposure from a frame with optional corrections.
cpl_error_code eris_ifu_fetch_std_param(const cpl_parameterlist *parlist, const char *recipename, struct stdParamStruct *stdParams)
Fetch standard parameters from parameter list into structure.
void eris_ifu_free_std_param(struct stdParamStruct *stdParams)
Free memory allocated for stdParamStruct.
cpl_vector * eris_ifu_image_collapse(const cpl_image *img)
cpl_error_code eris_ifu_calc_bpm(const cpl_parameterlist *pl, const char *recipe_name, hdrl_image *master_img, const hdrl_imagelist *imglist_on, cpl_mask **bpm2dMask, cpl_mask **bpm3dMask)
Create 2D and/or 3D bad pixel masks using HDRL algorithms.
cpl_error_code eris_ifu_add_badpix_border(cpl_image *data, cpl_boolean add_ones, cpl_image *dqi)
Flag detector border pixels as bad in image and optionally set to 1.
cpl_error_code eris_ifu_add_std_params(cpl_parameterlist *pl, const char *recipename)
Add standard recipe parameters to a parameter list.
cpl_error_code eris_ifu_exposure_column_correction(cpl_image *image)
Perform column correction on raw detector image.
cpl_vector * eris_ifu_polyfit_1d(const cpl_vector *x, const cpl_vector *y, const int degree)
Fit a 1D polynomial to vector data.
cpl_mask * eris_ifu_detect_crh(hdrl_image *image, int exposureCorrectionMode, hdrl_parameter *laCosmicParams, bool maskImage)
Detect cosmic ray hits using LA-cosmic algorithm.
cpl_error_code eris_parlist_config_add_bpm(cpl_parameterlist *pl, const char *recname)
Add bad pixel mask configuration parameters to parameter list.
cpl_image * eris_ifu_calc_noise_map(const cpl_image *data, double gain, double readnoise)
Calculate a noise map for detector data.
cpl_error_code eris_ifu_saturation_detection(cpl_image *image, cpl_image *dqi)
Detect saturated pixels and flag them as bad.
hdrl_image * eris_ifu_load_exposure_file(const char *filename, int exposureCorrectionMode, cpl_image *dqi)
Load a raw detector exposure from file with corrections and noise.
hdrl_imagelist * eris_ifu_load_exposure_frameset(const cpl_frameset *frameset, int exposureCorrectionMode)
Load raw exposures from a frameset with optional corrections.
ifsInstrument eris_ifu_get_instrument_frame(cpl_frame *frame)
Get instrument identifier from a CPL frame.
double eris_ifu_image_get_mean(const cpl_image *image)
‍**
hdrl_image * eris_ifu_raw_hdrl_image(const cpl_image *cplImage)
Create an HDRL image from a CPL image with calculated noise.
hdrl_image * eris_ifu_warp_polynomial_image(const hdrl_image *hdrlInImg, const cpl_polynomial *poly_u, const cpl_polynomial *poly_v)
Warp an HDRL image using 2D polynomial transformations.
void eris_ifu_free_propertylist(cpl_propertylist **item)
Free memory and set pointer to null.
void eris_ifu_free_string(char **item)
Free memory and set pointer to null.
cpl_error_code eris_ifu_parameterlist_append_list(cpl_parameterlist *p1, const cpl_parameterlist *p2)
Append a parameterlist to another one.
void eris_ifu_free_double_array(double **item)
Free memory and set pointer to null.
void eris_ifu_free_vector(cpl_vector **item)
Free memory and set pointer to null.
cpl_error_code eris_ifu_save_vector_dbg(const cpl_vector *vec, const char *filename, int create, const cpl_propertylist *pl)
Save vector for debugging (quick, no DFS overhead)
void eris_ifu_free_parameter(cpl_parameter **item)
Free memory and set pointer to null.
void eris_ifu_free_imagelist(cpl_imagelist **item)
Free memory and set pointer to null.
void eris_ifu_free_hdrl_image(hdrl_image **item)
Free memory and set pointer to null.
void eris_ifu_free_image(cpl_image **item)
Free memory and set pointer to null.
void eris_ifu_free_mask(cpl_mask **item)
Free memory and set pointer to null.
void eris_ifu_free_parameterlist(cpl_parameterlist **item)
Free memory and set pointer to null.
void eris_ifu_free_hdrl_parameter(hdrl_parameter **item)
Free memory and set pointer to null.
int eris_ifu_is_nan_or_inf(double A)
Checks if a value is nan, inf or -inf.
cpl_error_code eris_check_error_code(const char *func_id)
handle CPL errors
Definition: eris_utils.c:56
cpl_mask * hdrl_bpm_2d_compute(const hdrl_image *img_in, const hdrl_parameter *params)
Detect bad pixels on a single image with an iterative process.
Definition: hdrl_bpm_2d.c:1136
hdrl_parameter * hdrl_bpm_2d_parameter_parse_parlist(const cpl_parameterlist *parlist, const char *prefix)
Parse parameter list to create input parameters for the BPM_2D.
Definition: hdrl_bpm_2d.c:896
hdrl_parameter * hdrl_bpm_2d_parameter_create_filtersmooth(double kappa_low, double kappa_high, int maxiter, cpl_filter_mode filter, cpl_border_mode border, int smooth_x, int smooth_y)
Creates BPM_2D Parameters object for HDRL_BPM_2D_FILTERSMOOTH.
Definition: hdrl_bpm_2d.c:137
cpl_parameterlist * hdrl_bpm_2d_parameter_create_parlist(const char *base_context, const char *prefix, const char *method_def, const hdrl_parameter *filtersmooth_def, const hdrl_parameter *legendresmooth_def)
Create parameter list for the BPM_2D computation.
Definition: hdrl_bpm_2d.c:798
hdrl_parameter * hdrl_bpm_2d_parameter_create_legendresmooth(double kappa_low, double kappa_high, int maxiter, int steps_x, int steps_y, int filter_size_x, int filter_size_y, int order_x, int order_y)
Creates BPM_2D Parameters object for HDRL_BPM_2D_LEGENDRESMOOTH.
Definition: hdrl_bpm_2d.c:191
cpl_parameterlist * hdrl_bpm_3d_parameter_create_parlist(const char *base_context, const char *prefix, const hdrl_parameter *defaults)
Create a parameter list for the BPM_3D computation.
Definition: hdrl_bpm_3d.c:231
cpl_imagelist * hdrl_bpm_3d_compute(const hdrl_imagelist *imglist, const hdrl_parameter *params)
detect bad pixels on a stack of identical images
Definition: hdrl_bpm_3d.c:409
hdrl_parameter * hdrl_bpm_3d_parameter_parse_parlist(const cpl_parameterlist *parlist, const char *prefix)
Parse a parameterlist to create input parameters for the BPM_3D.
Definition: hdrl_bpm_3d.c:314
hdrl_parameter * hdrl_bpm_3d_parameter_create(double kappa_low, double kappa_high, hdrl_bpm_3d_method method)
Creates BPM Parameters object for the imagelist method.
Definition: hdrl_bpm_3d.c:108
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_mode_parameter_create(double histo_min, double histo_max, double bin_size, hdrl_mode_type mode_method, cpl_size error_niter)
create a parameter object for the mode
cpl_parameterlist * hdrl_collapse_parameter_create_parlist(const char *base_context, const char *prefix, const char *method_def, hdrl_parameter *sigclip_def, hdrl_parameter *minmax_def, hdrl_parameter *mode_def)
Create parameters for the collapse.
hdrl_parameter * hdrl_collapse_minmax_parameter_create(double nlow, double nhigh)
create a parameter object for min-max rejected mean
hdrl_parameter * hdrl_flat_parameter_create(cpl_size filter_size_x, cpl_size filter_size_y, hdrl_flat_method method)
Creates FLAT Parameters object.
Definition: hdrl_flat.c:173
cpl_parameterlist * hdrl_flat_parameter_create_parlist(const char *base_context, const char *prefix, const hdrl_parameter *defaults)
Create a parameter list for the FLAT computation.
Definition: hdrl_flat.c:258
cpl_error_code hdrl_image_reject_from_mask(hdrl_image *self, const cpl_mask *map)
set bpm of hdrl_image
Definition: hdrl_image.c:407
cpl_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
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
cpl_error_code hdrl_image_insert(hdrl_image *self, const cpl_image *image, const cpl_image *error, cpl_size xpos, cpl_size ypos)
Copy cpl images into an hdrl image.
Definition: hdrl_image.c:715
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
hdrl_image * hdrl_image_new(cpl_size nx, cpl_size ny)
create new zero filled hdrl image
Definition: hdrl_image.c:311
cpl_error_code hdrl_imagelist_set(hdrl_imagelist *himlist, hdrl_image *himg, cpl_size pos)
Insert an image into an imagelist.
cpl_size hdrl_imagelist_get_size(const hdrl_imagelist *himlist)
Get the number of images in the imagelist.
hdrl_imagelist * hdrl_imagelist_new(void)
Create an empty imagelist.
cpl_parameterlist * hdrl_lacosmic_parameter_create_parlist(const char *base_context, const char *prefix, const hdrl_parameter *defaults)
Create parameter list for the LaCosmic computation.
cpl_mask * hdrl_lacosmic_edgedetect(const hdrl_image *ima_in, const hdrl_parameter *params)
Detect bad-pixels / cosmic-rays on a single image.
hdrl_parameter * hdrl_lacosmic_parameter_parse_parlist(const cpl_parameterlist *parlist, const char *prefix)
Parse parameterlist to create input parameters for the LaCosmic.
hdrl_parameter * hdrl_lacosmic_parameter_create(double sigma_lim, double f_lim, int max_iter)
Creates LaCosmic parameters object.
void hdrl_parameter_delete(hdrl_parameter *obj)
shallow delete of a parameter