00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028 #ifdef HAVE_CONFIG_H
00029 #include <config.h>
00030 #endif
00031
00032
00033
00034
00035
00036 #ifdef HAVE_FFTW
00037 #include <complex.h>
00038 #include <fftw3.h>
00039 #endif
00040
00041 #include <math.h>
00042 #include <string.h>
00043 #include <assert.h>
00044 #include <float.h>
00045
00046 #include <cpl.h>
00047
00048 #include "irplib_detmon.h"
00049
00050 #include "irplib_hist.h"
00051 #include "irplib_utils.h"
00052
00053 #include "irplib_math_const.h"
00054
00055
00056 #define pdist(x1,y1,x2,y2) (((x1-x2)*(x1-x2))+((y1-y2)*(y1-y2)))
00057
00058 #define cpl_drand() ((double)rand()/(double)RAND_MAX)
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073 #define HIST_FACT 2.354820045
00074
00075 enum pixeltypes
00076 {
00077 HOT = 0,
00078 DEAD = 1,
00079 NOISY = 2
00080 };
00081
00082 enum stackingtypes
00083 {
00084 MINMAX = 0,
00085 MEAN = 1,
00086 MEDIAN = 2,
00087 KSIGMA = 3
00088 };
00089
00090 enum readouts
00091 {
00092 HORIZONTAL = 1,
00093 VERTICAL = 2
00094 };
00095
00096
00097 static struct
00098 {
00099
00100 const char *method;
00101 const char *pmethod;
00102 ronbias_method method_bitmask;
00103 int prescan_llx;
00104 int prescan_lly;
00105 int prescan_urx;
00106 int prescan_ury;
00107 int overscan_llx;
00108 int overscan_lly;
00109 int overscan_urx;
00110 int overscan_ury;
00111 int preoverscan_degree;
00112 int random_nsamples;
00113 int random_sizex;
00114 int random_sizey;
00115 int criteria;
00116 int ref_llx;
00117 int ref_lly;
00118 int ref_urx;
00119 int ref_ury;
00120 const char *stacking_method;
00121 int stacking_ks_low;
00122 int stacking_ks_high;
00123 int stacking_ks_iter;
00124 int master_shift_x;
00125 int master_shift_y;
00126 int ron_llx;
00127 int ron_lly;
00128 int ron_urx;
00129 int ron_ury;
00130 } detmon_ronbias_config;
00131
00132 static struct
00133 {
00134 int mode;
00135 cpl_boolean direction;
00136 double speed;
00137 int llx;
00138 int lly;
00139 int urx;
00140 int ury;
00141 int kappa;
00142 int exts;
00143 int nb_extensions;
00144 } detmon_pernoise_config;
00145
00146 static struct
00147 {
00148 const char * ron_method;
00149 const char * dsnu_method;
00150 int exts;
00151 int nb_extensions;
00152 cpl_boolean opt_nir;
00153 } detmon_dark_config;
00154
00155 #define NIR TRUE
00156 #define OPT FALSE
00157
00158
00159
00160
00161
00162
00163
00164
00165
00166 static cpl_error_code
00167 irplib_detmon_ronbias_retrieve_parlist(const char *,
00168 const char *,
00169 const cpl_parameterlist *);
00170
00171 static cpl_error_code
00172 irplib_detmon_ronbias_random(const cpl_imagelist *,
00173 const cpl_image *, cpl_propertylist *);
00174
00175 static cpl_error_code
00176 irplib_detmon_ronbias_histo(const cpl_imagelist *,
00177 const cpl_image *, cpl_propertylist *);
00178
00179 static cpl_error_code
00180 irplib_detmon_ronbias_preoverscan(const cpl_imagelist *,
00181 cpl_propertylist *, cpl_image **);
00182
00183 static cpl_error_code
00184 irplib_detmon_ronbias_region(const cpl_imagelist *,
00185 const cpl_image *, cpl_propertylist *);
00186
00187 static cpl_image *
00188 irplib_detmon_ronbias_master(const cpl_imagelist *,
00189 cpl_propertylist *);
00190
00191 static cpl_error_code
00192 irplib_detmon_ronbias_save(const cpl_parameterlist *,
00193 cpl_frameset *,
00194 const char *,
00195 const char *,
00196 const char *,
00197 const char *,
00198 const char *,
00199 const char *,
00200 const char *,
00201 const cpl_image *,
00202 const cpl_image *,
00203 cpl_propertylist *,
00204 const int, const int, cpl_frameset *);
00205
00206 int
00207 irplib_detmon_ronbias_dfs_set_groups(cpl_frameset *, const char *);
00208
00209 static cpl_image *
00210 irplib_detmon_build_synthetic(cpl_image *,
00211 cpl_image *);
00212
00213 static cpl_error_code
00214 irplib_detmon_ronbias_histo_reduce(const cpl_image * c_raw,
00215 double * bias,
00216 double * fwhm,
00217 double * max);
00218
00219 static cpl_error_code
00220 irplib_detmon_ronbias_dutycycl(const cpl_frameset *, cpl_propertylist *);
00221
00222 cpl_error_code
00223 irplib_detmon_rm_bpixs(cpl_image **,
00224 const double,
00225 int ,
00226 int );
00227
00228
00229
00230 static cpl_error_code
00231 irplib_flux_get_bias_window(const cpl_image *,
00232 const int *,
00233 int ,
00234 int ,
00235 double *,
00236 double *);
00237
00238 static cpl_bivector *
00239 irplib_bivector_gen_rect_poisson(const int *r,
00240 const int np,
00241 const int homog);
00242
00243
00244
00245
00246 cpl_error_code
00247 irplib_detmon_ronbias_check_defaults(const cpl_image *);
00248
00249
00250
00251
00252 int
00253 irplib_detmon_pernoise_dfs_set_groups(cpl_frameset *,
00254 const char *);
00255
00256 static cpl_error_code
00257 irplib_detmon_pernoise_retrieve_parlist(const char *,
00258 const char *,
00259 const cpl_parameterlist *);
00260
00261 static cpl_error_code
00262 irplib_detmon_pernoise_qc(cpl_propertylist *,
00263 cpl_table *,
00264 int);
00265
00266 static cpl_error_code
00267 irplib_detmon_pernoise_save(const cpl_parameterlist *,
00268 cpl_frameset *,
00269 const char *,
00270 const char *,
00271 const char *,
00272 const char *,
00273 cpl_table **,
00274 cpl_propertylist **,
00275 const int,
00276 const int,
00277 const cpl_frameset *);
00278
00279 cpl_error_code
00280 irplib_detmon_pernoise_rm_bg(cpl_image *,
00281 int,
00282 int);
00283
00284 int
00285 irplib_detmon_dark_dfs_set_groups(cpl_frameset *,
00286 const char *);
00287
00288 cpl_error_code
00289 irplib_detmon_dark_dsnu(cpl_frameset *,
00290 cpl_imagelist *,
00291 cpl_table *,
00292 cpl_image *,
00293 int pos);
00294
00295
00296 static cpl_error_code
00297 irplib_detmon_dark_save(const cpl_parameterlist *,
00298 cpl_frameset *,
00299 const char *,
00300 const char *,
00301 const char *,
00302 const char *,
00303 const char *,
00304 const char *,
00305 cpl_imagelist **,
00306 cpl_table **,
00307 cpl_imagelist **,
00308 cpl_propertylist **,
00309 const int,
00310 const int,
00311 const cpl_frameset *);
00312
00313
00314
00315 static cpl_error_code
00316 irplib_detmon_retrieve_dark_params(const char *,
00317 const char *,
00318 const cpl_parameterlist *);
00319
00320 cpl_error_code
00321 irplib_detmon_dark_qc(cpl_propertylist *,
00322 cpl_image *);
00323
00324
00325 #define CPL_NEW 0
00326 #define CPL_DELETE 1
00327
00328 #define CPL_TABLE 0
00329 #define CPL_PROPERTYLIST 1
00330 #define CPL_IMAGELIST 2
00331 #define CPL_IMAGE 3
00332
00333 #define CPL_OPERATION CPL_NEW
00334
00335 cpl_table ** irplib_detmon_table_new(int nexts)
00336 {
00337 #define CPL_CLASS CPL_TABLE
00338 #include "irplib_detmon_body.h"
00339 #undef CPL_CLASS
00340 }
00341
00342 cpl_propertylist ** irplib_detmon_propertylist_new(int nexts)
00343 {
00344 #define CPL_CLASS CPL_PROPERTYLIST
00345 #include "irplib_detmon_body.h"
00346 #undef CPL_CLASS
00347 }
00348
00349 cpl_imagelist ** irplib_detmon_imagelist_new(int nexts)
00350 {
00351 #define CPL_CLASS CPL_IMAGELIST
00352 #include "irplib_detmon_body.h"
00353 #undef CPL_CLASS
00354 }
00355
00356 cpl_image ** irplib_detmon_image_new(int nexts)
00357 {
00358 #define CPL_CLASS CPL_IMAGE
00359 #include "irplib_detmon_body.h"
00360 #undef CPL_CLASS
00361 }
00362
00363 #undef CPL_OPERATION
00364
00365 #define CPL_OPERATION CPL_DELETE
00366
00367 void irplib_detmon_table_delete(cpl_table ** d, int nexts)
00368 {
00369 #define CPL_CLASS CPL_TABLE
00370 #include "irplib_detmon_body.h"
00371 #undef CPL_CLASS
00372 }
00373
00374 void irplib_detmon_propertylist_delete(cpl_propertylist ** d, int nexts)
00375 {
00376 #define CPL_CLASS CPL_PROPERTYLIST
00377 #include "irplib_detmon_body.h"
00378 #undef CPL_CLASS
00379 }
00380
00381 void irplib_detmon_imagelist_delete(cpl_imagelist ** d, int nexts)
00382 {
00383 #define CPL_CLASS CPL_IMAGELIST
00384 #include "irplib_detmon_body.h"
00385 #undef CPL_CLASS
00386 }
00387
00388 void irplib_detmon_image_delete(cpl_image ** d, int nexts)
00389 {
00390 #define CPL_CLASS CPL_IMAGE
00391 #include "irplib_detmon_body.h"
00392 #undef CPL_CLASS
00393 }
00394
00395 #undef CPL_OPERATION
00396
00397 cpl_error_code
00398 irplib_detmon_fill_ronbias_params(cpl_parameterlist * parlist,
00399 const char *recipe_name,
00400 const char *pipeline_name)
00401 {
00402
00403 irplib_detmon_fill_parlist(parlist, recipe_name, pipeline_name, 29,
00404 "method",
00405 "Method to be used when computing bias. Methods appliable: <RANDOM | HISTO | PREOVERSCAN | REGION | ALL>. By default ALL methods are applied. More than a method can be chosen; in that case selected methods must be separated by a single space and put together between inverted commas (ex. --method=\"HISTO REGION\").\n RANDOM: Bias is computed as the mean value on a given number (--random.nsamples) of boxes (dimensions --random.sizex and --random.sizey) randomly taken accross the detector.\n HISTO: An histogram of the pixels of the image is built.\n PREOVERSCAN: Mean, median and RMS values computed and designated areas. \n REGION: Mean, median and RMS values on reference region.",
00406 "CPL_TYPE_STRING", "ALL",
00407 "pmethod",
00408 "Pre-method for RANDOM, HISTO and REGION."
00409 "Difference raw frames or not",
00410 "CPL_TYPE_STRING", "NORM",
00411 "prescan.llx",
00412 "x coordinate of the lower-left point of "
00413 "the prescan area",
00414 "CPL_TYPE_INT", -1,
00415 "prescan.lly",
00416 "y coordinate of the lower-left point of "
00417 "the prescan area",
00418 "CPL_TYPE_INT", -1,
00419 "prescan.urx",
00420 "x coordinate of the upper-right point of "
00421 "the prescan area",
00422 "CPL_TYPE_INT", -1,
00423 "prescan.ury",
00424 "y coordinate of the upper-right point of "
00425 "the prescan area",
00426 "CPL_TYPE_INT", -1,
00427 "overscan.llx",
00428 "x coordinate of the lower-left point of "
00429 "the prescan area",
00430 "CPL_TYPE_INT", -1,
00431 "overscan.lly",
00432 "y coordinate of the lower-left point of "
00433 "the overscan area",
00434 "CPL_TYPE_INT", -1,
00435 "overscan.urx",
00436 "x coordinate of the upper-right point of "
00437 "the overscan area",
00438 "CPL_TYPE_INT", -1,
00439 "overscan.ury",
00440 "y coordinate of the upper-right point of "
00441 "the overscan area",
00442 "CPL_TYPE_INT", -1,
00443 "preoverscan.degree",
00444 "Degree used for pre-overscan method",
00445 "CPL_TYPE_INT", 1,
00446 "random.nsamples",
00447 "Number of samples",
00448 "CPL_TYPE_INT", 50,
00449 "random.sizex",
00450 "X size of the boxes",
00451 "CPL_TYPE_INT", 100,
00452 "random.sizey",
00453 "Y size of the boxes",
00454 "CPL_TYPE_INT", 100,
00455 "criteria",
00456 "Criteria",
00457 "CPL_TYPE_INT", 0,
00458 "ref.llx",
00459 "x coordinate of the lower-left point "
00460 "of the reference region of the frame",
00461 "CPL_TYPE_INT", -1,
00462 "ref.lly",
00463 "y coordinate of the lower-left point "
00464 "of the reference region of the frame",
00465 "CPL_TYPE_INT", -1,
00466 "ref.urx",
00467 "x coordinate of the upper-right point "
00468 "of the reference region of the frame",
00469 "CPL_TYPE_INT", -1,
00470 "ref.ury",
00471 "y coordinate of the upper-right point "
00472 "of the reference region of the frame",
00473 "CPL_TYPE_INT", -1,
00474 "stacking.method",
00475 "Method to be used when stacking the master. Posible values < MINMAX | MEAN | MEDIAN | KSIGMA >",
00476 "CPL_TYPE_STRING", "MEAN",
00477 "stacking.ks.low",
00478 "Low threshold for kappa-sigma clipping",
00479 "CPL_TYPE_INT", 3,
00480 "stacking.ks.high",
00481 "High threshold for kappa-sigma clipping",
00482 "CPL_TYPE_INT", 10000,
00483 "stacking.ks.iter",
00484 "Nb of iterations for kappa-sigma clipping",
00485 "CPL_TYPE_INT", 25,
00486 "master.shift.x",
00487 "Master shift X",
00488 "CPL_TYPE_INT", 0,
00489 "master.shift.y",
00490 "Master shift Y",
00491 "CPL_TYPE_INT", 0,
00492 "ron.llx",
00493 "x coordinate of the lower-left point "
00494 "of the RON frame",
00495 "CPL_TYPE_INT", -1,
00496 "ron.lly",
00497 "y coordinate of the lower-left point "
00498 "of the RON frame",
00499 "CPL_TYPE_INT", -1,
00500 "ron.urx",
00501 "x coordinate of the upper-right point "
00502 "of the RON frame",
00503 "CPL_TYPE_INT", -1,
00504 "ron.ury",
00505 "y coordinate of the upper-right point "
00506 "of the RON frame", "CPL_TYPE_INT", -1);
00507
00508 return 0;
00509 }
00510
00511 int
00512 irplib_detmon_fill_parlist(cpl_parameterlist * parlist,
00513 const char *recipe_name, const char *pipeline_name,
00514 int npars, ...)
00515 {
00516
00517 va_list ap;
00518
00519 char *group_name;
00520
00521 int pars_counter = 0;
00522
00523 group_name = cpl_sprintf("%s.%s", pipeline_name, recipe_name);
00524 assert(group_name != NULL);
00525
00526 #define insert_par(PARNAME, PARDESC, PARVALUE, PARTYPE) \
00527 do { \
00528 char * par_name = cpl_sprintf("%s.%s", group_name, PARNAME); \
00529 cpl_parameter * p; \
00530 assert(par_name != NULL); \
00531 p = cpl_parameter_new_value(par_name, PARTYPE, \
00532 PARDESC, group_name, PARVALUE); \
00533 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, PARNAME); \
00534 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV); \
00535 cpl_parameterlist_append(parlist, p); \
00536 cpl_free(par_name); \
00537 } while(0);
00538
00539
00540 va_start(ap, npars);
00541
00542 while(pars_counter < npars) {
00543 char *name = va_arg(ap, char *);
00544 char *desc = va_arg(ap, char *);
00545 char *type = va_arg(ap, char *);
00546
00547 if(!strcmp(type, "CPL_TYPE_INT")) {
00548 int v1 = va_arg(ap, int);
00549
00550 insert_par(name, desc, v1, CPL_TYPE_INT);
00551 } else if(!strcmp(type, "CPL_TYPE_BOOL")) {
00552 char *v2 = va_arg(ap, char *);
00553
00554 if(!strcmp(v2, "CPL_FALSE"))
00555 insert_par(name, desc, CPL_FALSE, CPL_TYPE_BOOL);
00556 if(!strcmp(v2, "CPL_TRUE"))
00557 insert_par(name, desc, CPL_TRUE, CPL_TYPE_BOOL);
00558 } else if(!strcmp(type, "CPL_TYPE_STRING")) {
00559 char *v2 = va_arg(ap, char *);
00560
00561 insert_par(name, desc, v2, CPL_TYPE_STRING);
00562 } else if(!strcmp(type, "CPL_TYPE_DOUBLE")) {
00563 double v3 = va_arg(ap, double);
00564 insert_par(name, desc, v3, CPL_TYPE_DOUBLE);
00565 }
00566
00567 pars_counter++;
00568 }
00569
00570 va_end(ap);
00571
00572 cpl_free(group_name);
00573
00574 #undef insert_par
00575 return 0;
00576 }
00577
00578
00579
00580
00581
00582
00583
00584
00585
00586
00587
00588
00589 int
00590 irplib_detmon_retrieve_par(const char *parn,
00591 const char *pipeline_name,
00592 const char *recipe_name,
00593 const cpl_parameterlist * parlist)
00594 {
00595 char *par_name;
00596 cpl_parameter *par;
00597 int value;
00598
00599 par_name = cpl_sprintf("%s.%s.%s", pipeline_name, recipe_name, parn);
00600 assert(par_name != NULL);
00601 par = cpl_parameterlist_find((cpl_parameterlist *) parlist, par_name);
00602 value = cpl_parameter_get_int(par);
00603 cpl_free(par_name);
00604
00605 return value;
00606 }
00607
00608
00609
00610
00611
00612
00613
00614
00615
00616
00617
00618
00619
00620 int
00621 irplib_detmon_compare_dits(const cpl_frame * frame1, const cpl_frame * frame2)
00622 {
00623 int comparison;
00624 cpl_propertylist *plist1;
00625 cpl_propertylist *plist2;
00626 double dval1, dval2;
00627
00628
00629 if(frame1 == NULL || frame2 == NULL)
00630 return -1;
00631
00632
00633 if((plist1 = cpl_propertylist_load(cpl_frame_get_filename(frame1),
00634 0)) == NULL) {
00635 cpl_msg_error(cpl_func, "getting header from reference frame");
00636 return -1;
00637 }
00638 if((plist2 = cpl_propertylist_load(cpl_frame_get_filename(frame2),
00639 0)) == NULL) {
00640 cpl_msg_error(cpl_func, "getting header from reference frame");
00641 cpl_propertylist_delete(plist1);
00642 return -1;
00643 }
00644
00645
00646 if(cpl_error_get_code()) {
00647 cpl_propertylist_delete(plist1);
00648 cpl_propertylist_delete(plist2);
00649 return -1;
00650 }
00651
00652
00653 comparison = 1;
00654 dval1 = irplib_pfits_get_exptime(plist1);
00655 dval2 = irplib_pfits_get_exptime(plist2);
00656 if(cpl_error_get_code()) {
00657 cpl_msg_error(cpl_func, "cannot get exposure time");
00658 cpl_propertylist_delete(plist1);
00659 cpl_propertylist_delete(plist2);
00660 return -1;
00661 }
00662 if(fabs(dval1 - dval2) > 1e-3)
00663 comparison = 0;
00664
00665
00666 cpl_propertylist_delete(plist1);
00667 cpl_propertylist_delete(plist2);
00668 return comparison;
00669 }
00670
00671
00672 double
00673 irplib_pfits_get_exptime(const cpl_propertylist * plist)
00674 {
00675 double exptime;
00676
00677 exptime = cpl_propertylist_get_double(plist, "EXPTIME");
00678
00679 return exptime;
00680 }
00681
00682
00683
00684
00685
00686
00687
00688
00689
00690
00691
00692
00693 static cpl_error_code
00694 irplib_detmon_ronbias_retrieve_parlist(const char *pipeline_name,
00695 const char *recipe_name,
00696 const cpl_parameterlist * parlist)
00697 {
00698 char *par_name;
00699 cpl_parameter *par;
00700
00701 char m1[20] = "";
00702 char m2[20] = "";
00703 char m3[20] = "";
00704
00705
00706 par_name = cpl_sprintf("%s.%s.method", pipeline_name, recipe_name);
00707 assert(par_name != NULL);
00708 par = cpl_parameterlist_find((cpl_parameterlist *) parlist, par_name);
00709 detmon_ronbias_config.method = cpl_parameter_get_string(par);
00710 cpl_free(par_name);
00711
00712 detmon_ronbias_config.method_bitmask = 0;
00713
00714 sscanf(detmon_ronbias_config.method, "%s %s %s", m1, m2, m3);
00715
00716 if(!strcmp(m1, "RANDOM") || !strcmp(m2, "RANDOM")
00717 || !strcmp(m3, "RANDOM"))
00718 detmon_ronbias_config.method_bitmask += RANDOM;
00719
00720 if(!strcmp(m1, "HISTO") || !strcmp(m2, "HISTO") || !strcmp(m3, "HISTO"))
00721 detmon_ronbias_config.method_bitmask += HISTO;
00722
00723 if(!strcmp(m1, "PREOVERSCAN") || !strcmp(m2, "PREOVERSCAN")
00724 || !strcmp(m3, "PREOVERSCAN"))
00725 detmon_ronbias_config.method_bitmask += PREOVERSCAN;
00726
00727 if(!strcmp(m1, "REGION") || !strcmp(m2, "REGION")
00728 || !strcmp(m3, "REGION"))
00729 detmon_ronbias_config.method_bitmask += REGION;
00730
00731 if(!strcmp(m1, "ALL"))
00732 detmon_ronbias_config.method_bitmask =
00733 RANDOM | HISTO | PREOVERSCAN | REGION;
00734
00735
00736 par_name = cpl_sprintf("%s.%s.pmethod", pipeline_name, recipe_name);
00737 assert(par_name != NULL);
00738 par = cpl_parameterlist_find((cpl_parameterlist *) parlist, par_name);
00739 detmon_ronbias_config.pmethod = cpl_parameter_get_string(par);
00740 cpl_free(par_name);
00741
00742
00743 detmon_ronbias_config.prescan_llx =
00744 irplib_detmon_retrieve_par("prescan.llx", pipeline_name, recipe_name,
00745 parlist);
00746
00747
00748 detmon_ronbias_config.prescan_lly =
00749 irplib_detmon_retrieve_par("prescan.lly", pipeline_name, recipe_name,
00750 parlist);
00751
00752 detmon_ronbias_config.prescan_urx =
00753 irplib_detmon_retrieve_par("prescan.urx", pipeline_name, recipe_name,
00754 parlist);
00755
00756 detmon_ronbias_config.prescan_ury =
00757 irplib_detmon_retrieve_par("prescan.ury", pipeline_name, recipe_name,
00758 parlist);
00759
00760 detmon_ronbias_config.overscan_llx =
00761 irplib_detmon_retrieve_par("overscan.llx", pipeline_name, recipe_name,
00762 parlist);
00763
00764 detmon_ronbias_config.overscan_lly =
00765 irplib_detmon_retrieve_par("overscan.lly", pipeline_name, recipe_name,
00766 parlist);
00767
00768 detmon_ronbias_config.overscan_urx =
00769 irplib_detmon_retrieve_par("overscan.urx", pipeline_name, recipe_name,
00770 parlist);
00771
00772 detmon_ronbias_config.overscan_ury =
00773 irplib_detmon_retrieve_par("overscan.ury", pipeline_name, recipe_name,
00774 parlist);
00775
00776
00777 detmon_ronbias_config.preoverscan_degree =
00778 irplib_detmon_retrieve_par("preoverscan.degree", pipeline_name,
00779 recipe_name, parlist);
00780
00781
00782 detmon_ronbias_config.random_nsamples =
00783 irplib_detmon_retrieve_par("random.nsamples", pipeline_name,
00784 recipe_name, parlist);
00785
00786
00787 detmon_ronbias_config.random_sizex =
00788 irplib_detmon_retrieve_par("random.sizex", pipeline_name,
00789 recipe_name, parlist);
00790
00791
00792 detmon_ronbias_config.random_sizey =
00793 irplib_detmon_retrieve_par("random.sizey", pipeline_name,
00794 recipe_name, parlist);
00795
00796
00797 detmon_ronbias_config.criteria =
00798 irplib_detmon_retrieve_par("criteria", pipeline_name, recipe_name,
00799 parlist);
00800
00801
00802 detmon_ronbias_config.ref_llx =
00803 irplib_detmon_retrieve_par("ref.llx", pipeline_name, recipe_name,
00804 parlist);
00805
00806 detmon_ronbias_config.ref_lly =
00807 irplib_detmon_retrieve_par("ref.lly", pipeline_name, recipe_name,
00808 parlist);
00809
00810 detmon_ronbias_config.ref_urx =
00811 irplib_detmon_retrieve_par("ref.urx", pipeline_name, recipe_name,
00812 parlist);
00813
00814 detmon_ronbias_config.ref_ury =
00815 irplib_detmon_retrieve_par("ref.ury", pipeline_name, recipe_name,
00816 parlist);
00817
00818
00819 par_name =
00820 cpl_sprintf("%s.%s.stacking.method", pipeline_name, recipe_name);
00821 assert(par_name != NULL);
00822 par = cpl_parameterlist_find((cpl_parameterlist *) parlist, par_name);
00823 detmon_ronbias_config.stacking_method = cpl_parameter_get_string(par);
00824 cpl_free(par_name);
00825
00826
00827 detmon_ronbias_config.stacking_ks_low =
00828 irplib_detmon_retrieve_par("stacking.ks.low", pipeline_name,
00829 recipe_name, parlist);
00830
00831 detmon_ronbias_config.stacking_ks_high =
00832 irplib_detmon_retrieve_par("stacking.ks.high", pipeline_name,
00833 recipe_name, parlist);
00834
00835 detmon_ronbias_config.stacking_ks_iter =
00836 irplib_detmon_retrieve_par("stacking.ks.iter", pipeline_name,
00837 recipe_name, parlist);
00838
00839 detmon_ronbias_config.master_shift_x =
00840 irplib_detmon_retrieve_par("master.shift.x", pipeline_name,
00841 recipe_name, parlist);
00842
00843 detmon_ronbias_config.master_shift_y =
00844 irplib_detmon_retrieve_par("master.shift.y", pipeline_name,
00845 recipe_name, parlist);
00846
00847 detmon_ronbias_config.ron_llx =
00848 irplib_detmon_retrieve_par("ron.llx", pipeline_name, recipe_name,
00849 parlist);
00850
00851 detmon_ronbias_config.ron_lly =
00852 irplib_detmon_retrieve_par("ron.lly", pipeline_name, recipe_name,
00853 parlist);
00854
00855 detmon_ronbias_config.ron_urx =
00856 irplib_detmon_retrieve_par("ron.urx", pipeline_name, recipe_name,
00857 parlist);
00858
00859 detmon_ronbias_config.ron_ury =
00860 irplib_detmon_retrieve_par("ron.ury", pipeline_name, recipe_name,
00861 parlist);
00862
00863 if(cpl_error_get_code()) {
00864 cpl_msg_error(cpl_func, "Failed to retrieve the input parameters");
00865 cpl_ensure_code(0, CPL_ERROR_DATA_NOT_FOUND);
00866 }
00867
00868
00869 return CPL_ERROR_NONE;
00870 }
00871
00872
00873
00874
00875
00876
00877
00878
00879
00880
00881
00882
00883 cpl_error_code
00884 irplib_detmon_ronbias_check_defaults(const cpl_image * reference)
00885 {
00886 const int nx = cpl_image_get_size_x(reference);
00887 const int ny = cpl_image_get_size_y(reference);
00888 cpl_error_code error;
00889
00890 if(detmon_ronbias_config.prescan_llx == -1)
00891 detmon_ronbias_config.prescan_llx = 1;
00892 if(detmon_ronbias_config.prescan_lly == -1)
00893 detmon_ronbias_config.prescan_lly = 1;
00894 if(detmon_ronbias_config.prescan_urx == -1)
00895 detmon_ronbias_config.prescan_urx = nx;
00896 if(detmon_ronbias_config.prescan_ury == -1)
00897 detmon_ronbias_config.prescan_ury = ny / 8;
00898
00899 if(detmon_ronbias_config.overscan_llx == -1)
00900 detmon_ronbias_config.overscan_llx = nx * 7 / 8;
00901 if(detmon_ronbias_config.overscan_lly == -1)
00902 detmon_ronbias_config.overscan_lly = 1;
00903 if(detmon_ronbias_config.overscan_urx == -1)
00904 detmon_ronbias_config.overscan_urx = nx;
00905 if(detmon_ronbias_config.overscan_ury == -1)
00906 detmon_ronbias_config.overscan_ury = ny;
00907
00908 if(detmon_ronbias_config.ref_llx == -1)
00909 detmon_ronbias_config.ref_llx = nx / 8;
00910 if(detmon_ronbias_config.ref_lly == -1)
00911 detmon_ronbias_config.ref_lly = ny / 8;
00912 if(detmon_ronbias_config.ref_urx == -1)
00913 detmon_ronbias_config.ref_urx = nx * 7 / 8;
00914 if(detmon_ronbias_config.ref_ury == -1)
00915 detmon_ronbias_config.ref_ury = ny * 7 / 8;
00916
00917 if(detmon_ronbias_config.ron_llx == -1)
00918 detmon_ronbias_config.ron_llx = 1;
00919 if(detmon_ronbias_config.ron_lly == -1)
00920 detmon_ronbias_config.ron_lly = 1;
00921 if(detmon_ronbias_config.ron_urx == -1)
00922 detmon_ronbias_config.ron_urx = nx;
00923 if(detmon_ronbias_config.ron_ury == -1)
00924 detmon_ronbias_config.ron_ury = ny;
00925
00926 error = cpl_error_get_code();
00927 cpl_ensure_code(!error, error);
00928
00929 return CPL_ERROR_NONE;
00930 }
00931
00946 cpl_error_code
00947 irplib_get_clean_mean_window(const cpl_image * img,
00948 const int llx,
00949 const int lly,
00950 const int urx, int ury,
00951 const int kappa,
00952 const int nclip,
00953 double *clean_mean, double *clean_stdev)
00954 {
00955
00956
00957 double mean = 0;
00958 double stdev = 0;
00959 double threshold = 0;
00960 double lo_cut = 0;
00961 double hi_cut = 0;
00962 double lo_cut_p = 0;
00963 double hi_cut_p = 0;
00964
00965 cpl_mask *mask = NULL;
00966 cpl_image *tmp = NULL;
00967 cpl_stats *stats = NULL;
00968 int i = 0;
00969 cpl_error_code error;
00970
00971 tmp = cpl_image_extract(img, llx, lly, urx, ury);
00972 cpl_ensure_code(tmp != NULL, cpl_error_get_code());
00973
00974 error = cpl_image_accept_all(tmp);
00975 cpl_ensure_code(!error, error);
00976
00977 for(i = 0; i < nclip; i++) {
00978
00979
00980 cpl_stats_delete(stats);
00981 stats =
00982 cpl_stats_new_from_image(tmp, CPL_STATS_MEAN | CPL_STATS_STDEV);
00983 cpl_ensure_code(stats != NULL, cpl_error_get_code());
00984 mean = cpl_stats_get_mean(stats);
00985 stdev = cpl_stats_get_stdev(stats);
00986
00987 threshold = kappa * stdev;
00988 lo_cut = mean - threshold;
00989 hi_cut = mean + threshold;
00990
00991 error = cpl_image_accept_all(tmp);
00992 cpl_ensure_code(!error, error);
00993 mask = cpl_mask_threshold_image_create(tmp, lo_cut, hi_cut);
00994 cpl_ensure_code(mask != NULL, cpl_error_get_code());
00995
00996 error = cpl_mask_not(mask);
00997 cpl_ensure_code(!error, error);
00998
00999 error = cpl_image_reject_from_mask(tmp, mask);
01000 cpl_ensure_code(!error, error);
01001 cpl_mask_delete(mask);
01002
01003 if(lo_cut_p == lo_cut && hi_cut_p == hi_cut)
01004 break;
01005 else {
01006 lo_cut_p = lo_cut;
01007 hi_cut_p = hi_cut;
01008 }
01009
01010 }
01011 *clean_mean = mean;
01012 *clean_stdev = stdev;
01013 cpl_image_delete(tmp);
01014 cpl_stats_delete(stats);
01015
01016 return 0;
01017
01018
01019 }
01020
01021
01022
01023
01024
01025
01026
01027
01028
01029
01030
01031
01032 cpl_error_code
01033 irplib_detmon_ronbias(cpl_frameset * frameset,
01034 const cpl_parameterlist * parlist,
01035 const char *tag,
01036 const char *recipe_name,
01037 const char *pipeline_name,
01038 const char *procatg_master,
01039 const char *procatg_xstr,
01040 const char *procatg_ystr,
01041 const char *procatg_bpm,
01042 const char *package,
01043 int (*compare) (const cpl_frame *, const cpl_frame *))
01044 {
01045
01046 int nsets;
01047 int *selection = NULL;
01048 int i;
01049 cpl_error_code error;
01050
01051 if(irplib_detmon_ronbias_dfs_set_groups(frameset, tag)) {
01052 cpl_msg_error(cpl_func, "Cannot identify RAW and CALIB frames");
01053 }
01054
01055
01056
01057
01058
01059
01060
01061
01062
01063
01064
01065
01066
01067
01068 irplib_detmon_ronbias_retrieve_parlist(pipeline_name,
01069 recipe_name, parlist);
01070
01071
01072
01073 if(compare == NULL)
01074 nsets = 1;
01075 else {
01076 cpl_msg_info(cpl_func, "Identify the different settings");
01077 selection = cpl_frameset_labelise(frameset, compare, &nsets);
01078 if(selection == NULL)
01079 cpl_msg_error(cpl_func, "Cannot labelise input frames");
01080 }
01081
01082
01083 for(i = 0; i < nsets; i++) {
01084 cpl_propertylist *qclist = cpl_propertylist_new();
01085 cpl_image *synthetic = NULL;
01086 cpl_image *masterbias;
01087
01088 cpl_imagelist *rawbiases
01089 = cpl_imagelist_load_frameset(frameset, CPL_TYPE_FLOAT, 1, 0);
01090 cpl_ensure_code(rawbiases != NULL, cpl_error_get_code());
01091
01092 error =
01093 irplib_detmon_ronbias_check_defaults(cpl_imagelist_get_const
01094 (rawbiases, 0));
01095 cpl_ensure_code(!error, error);
01096
01097 error = irplib_detmon_ronbias_dutycycl(frameset, qclist);
01098 cpl_ensure_code(!error, error);
01099
01100 masterbias = irplib_detmon_ronbias_master(rawbiases, qclist);
01101
01102
01103
01104
01105
01106
01107
01108
01109
01110 cpl_ensure_code(masterbias != NULL, CPL_ERROR_DATA_NOT_FOUND);
01111 cpl_ensure_code(!cpl_error_get_code(), cpl_error_get_code());
01112
01113
01114
01115
01116
01117
01118 if(detmon_ronbias_config.method_bitmask & RANDOM) {
01119
01120 irplib_detmon_ronbias_random(rawbiases, masterbias, qclist);
01121 cpl_ensure_code(!cpl_error_get_code(), cpl_error_get_code());
01122 }
01123
01124 if(detmon_ronbias_config.method_bitmask & HISTO) {
01125 irplib_detmon_ronbias_histo(rawbiases, masterbias, qclist);
01126 cpl_ensure_code(!cpl_error_get_code(), cpl_error_get_code());
01127 }
01128
01129 if(detmon_ronbias_config.method_bitmask & PREOVERSCAN) {
01130 irplib_detmon_ronbias_preoverscan(rawbiases,
01131 qclist, &synthetic);
01132 cpl_ensure_code(!cpl_error_get_code(), cpl_error_get_code());
01133 }
01134
01135 if(detmon_ronbias_config.method_bitmask & REGION) {
01136 irplib_detmon_ronbias_region(rawbiases, masterbias,
01137 qclist);
01138 cpl_ensure_code(!cpl_error_get_code(), cpl_error_get_code());
01139 }
01140
01141
01142
01143
01144
01145
01146
01147
01148 irplib_detmon_ronbias_save(parlist, frameset, recipe_name,
01149 pipeline_name, procatg_master,
01150 procatg_xstr, procatg_ystr, procatg_bpm,
01151 package, masterbias, synthetic, qclist, 0,
01152 0, frameset);
01153
01154 cpl_propertylist_delete(qclist);
01155 cpl_image_delete(synthetic);
01156 cpl_image_delete(masterbias);
01157 cpl_imagelist_delete(rawbiases);
01158 }
01159
01160 return cpl_error_get_code();
01161 }
01162
01163
01164
01165
01166
01167
01168
01169
01170
01171
01172
01173
01174 static cpl_error_code
01175 irplib_detmon_ronbias_random(const cpl_imagelist * rawbiases,
01176 const cpl_image * masterbias,
01177 cpl_propertylist * qclist)
01178 {
01179 int nraws = cpl_imagelist_get_size(rawbiases);
01180 int i;
01181 double bias = DBL_MAX;
01182 double bias_error;
01183
01184 double ron_error;
01185 double *ron =
01186 (double *) cpl_malloc(sizeof(double) * nraws);
01187
01188 cpl_vector *v;
01189 double stdev;
01190 cpl_error_code error;
01191
01192
01193
01194 if(!strcmp(detmon_ronbias_config.pmethod, "DIF")) {
01195 nraws--;
01196
01197
01198 for(i = 0; i < nraws; i++) {
01199 const cpl_image *c1_raw =
01200 cpl_imagelist_get_const(rawbiases, i);
01201 const cpl_image *c2_raw =
01202 cpl_imagelist_get_const(rawbiases, i + 1);
01203
01204 const cpl_image *c_raw = cpl_image_subtract_create(c1_raw,
01205 c2_raw);
01206 cpl_flux_get_noise_window(c_raw, NULL,
01207 detmon_ronbias_config.random_sizex / 2,
01208 detmon_ronbias_config.random_nsamples,
01209 ron + i, &ron_error);
01210 cpl_image_delete((cpl_image*)c_raw);
01211
01212 }
01213 } else {
01214 for(i = 0; i < nraws; i++) {
01215 const cpl_image *c_raw = cpl_imagelist_get_const(rawbiases, i);
01216
01217 cpl_flux_get_noise_window(c_raw, NULL,
01218 detmon_ronbias_config.random_sizex / 2,
01219 detmon_ronbias_config.random_nsamples,
01220 ron + i, &ron_error);
01221 }
01222 }
01223
01224
01225
01226
01227 irplib_flux_get_bias_window(masterbias, NULL,
01228 detmon_ronbias_config.random_sizex / 2,
01229 detmon_ronbias_config.random_nsamples,
01230 &bias, &bias_error);
01231
01232 v = cpl_vector_wrap(nraws, ron);
01233 stdev = cpl_vector_get_median_const(v);
01234 cpl_vector_unwrap(v);
01235
01236 cpl_free(ron);
01237
01238 error = cpl_propertylist_append_double(qclist, "ESO QC BIAS RANDOM BIAS", bias);
01239 cpl_ensure_code(!error, error);
01240 error = cpl_propertylist_append_double(qclist, "ESO QC BIAS RANDOM RMS", stdev);
01241 cpl_ensure_code(!error, error);
01242 return cpl_error_get_code();
01243 }
01244
01245 static cpl_error_code
01246 irplib_detmon_ronbias_histo(const cpl_imagelist * rawbiases,
01247 const cpl_image * masterbias,
01248 cpl_propertylist * qclist)
01249 {
01250 int nraws = cpl_imagelist_get_size(rawbiases);
01251 int i;
01252
01253 double mbias, mfwhm, mmax;
01254
01255 cpl_vector * fwhms;
01256 cpl_vector * maxs;
01257 cpl_error_code error;
01258
01259 if(!strcmp(detmon_ronbias_config.pmethod, "DIF")) nraws--;
01260
01261 fwhms = cpl_vector_new(nraws);
01262 maxs = cpl_vector_new(nraws);
01263
01264 for(i = 0; i < nraws; i++) {
01265
01266 const cpl_image * c_raw;
01267 double bias, fwhm, max;
01268
01269 if(strcmp(detmon_ronbias_config.pmethod, "DIF")) {
01270 c_raw = cpl_imagelist_get_const(rawbiases, i);
01271 } else {
01272 const cpl_image *c1_raw = cpl_imagelist_get_const(rawbiases, i);
01273 const cpl_image * c2_raw = cpl_imagelist_get_const(rawbiases, i+1);
01274 c_raw = cpl_image_subtract_create(c1_raw, c2_raw);
01275 }
01276
01277 error = irplib_detmon_ronbias_histo_reduce(c_raw, &bias, &fwhm, &max);
01278 cpl_ensure_code(!error, error);
01279
01280 if(!strcmp(detmon_ronbias_config.pmethod, "DIF"))
01281 cpl_image_delete((cpl_image *)c_raw);
01282
01283 error = cpl_vector_set(maxs, i, max);
01284 cpl_ensure_code(!error, error);
01285
01286 error = cpl_vector_set(fwhms, i, fwhm);
01287 cpl_ensure_code(!error, error);
01288
01289
01290
01291
01292
01293
01294 }
01295
01296 error = cpl_vector_divide_scalar(fwhms, HIST_FACT);
01297 cpl_ensure_code(!error, error);
01298
01299 irplib_detmon_ronbias_histo_reduce(masterbias, &mbias, &mfwhm, &mmax);
01300
01301 error = cpl_propertylist_append_double(qclist, "ESO QC BIAS HISTO BIAS",
01302 mbias);
01303 cpl_ensure_code(!error, error);
01304 error = cpl_propertylist_append_double(qclist, "ESO QC BIAS HISTO RMS",
01305 cpl_vector_get_mean(fwhms));
01306 cpl_ensure_code(!error, error);
01307 error =
01308 cpl_propertylist_append_double(qclist, "ESO QC BIAS HISTO RMS STDEV",
01309 cpl_vector_get_stdev(fwhms));
01310 cpl_ensure_code(!error, error);
01311 error = cpl_propertylist_append_double(qclist, "ESO QC BIAS HISTO MAX",
01312 cpl_vector_get_mean(maxs));
01313 cpl_ensure_code(!error, error);
01314 error =
01315 cpl_propertylist_append_double(qclist, "ESO QC BIAS HISTO MAX STDEV",
01316 cpl_vector_get_stdev(maxs));
01317 cpl_ensure_code(!error, error);
01318
01319 cpl_vector_delete(fwhms);
01320 cpl_vector_delete(maxs);
01321
01322 return CPL_ERROR_NONE;
01323 }
01324
01325 static cpl_error_code
01326 irplib_detmon_ronbias_histo_reduce(const cpl_image * c_raw,
01327 double * bias,
01328 double * fwhm,
01329 double * max)
01330 {
01331 int j, k;
01332 unsigned long uj;
01333 int rej;
01334 irplib_hist *hist;
01335 unsigned long maxwhere = 0;
01336 unsigned long maxf;
01337
01338 double mean, stdev;
01339 cpl_image * dupi;
01340 unsigned long x1a = 1;
01341 unsigned long x2a = 1;
01342
01343 double x1 = 0;
01344 double x2 = 0;
01345
01346 double maxwhere_interp;
01347 double max_interp;
01348 double a, b, c;
01349 cpl_matrix * coeffs =cpl_matrix_new(3, 3);
01350 cpl_matrix * rhs =cpl_matrix_new(3, 1);
01351 int p, q;
01352 cpl_matrix * result = NULL;
01353
01354
01355 int points = 0;
01356
01357 cpl_error_code error;
01358 mean = cpl_image_get_mean(c_raw);
01359 stdev = cpl_image_get_stdev(c_raw);
01360 dupi = cpl_image_duplicate(c_raw);
01361
01362
01363
01364
01365
01366
01367 for(j = 1; j <= cpl_image_get_size_x(c_raw); j++) {
01368 for(k = 1; k <= cpl_image_get_size_y(c_raw); k++) {
01369 if(fabs(cpl_image_get(c_raw, j, k, &rej) - mean) > stdev * 3)
01370 cpl_image_set(dupi, j, k, mean);
01371 }
01372 }
01373
01374
01375 hist = irplib_hist_new();
01376 error = irplib_hist_fill(hist, dupi);
01377 cpl_ensure_code(!error, error);
01378
01379 cpl_image_delete(dupi);
01380
01381 maxf = irplib_hist_get_max(hist, &maxwhere);
01382
01383 for( p = 0; p< 3; p++){
01384 unsigned long bi = irplib_hist_get_value(hist, maxwhere-1+p);
01385 cpl_matrix_set(rhs, p, 0, bi);
01386 for( q= 0; q< 3; q++) {
01387 cpl_matrix_set(coeffs, p,q,pow((maxwhere-1+p),q));
01388 }
01389 }
01390
01391 result = cpl_matrix_solve(coeffs, rhs);
01392
01393 a = cpl_matrix_get(result, 2, 0);
01394 b = cpl_matrix_get(result, 1, 0);
01395 c = cpl_matrix_get(result, 0, 0);
01396
01397 maxwhere_interp = -0.5 * b / (2 * a);
01398 max_interp = -1 * b * b / (4 * a) + c;
01399
01400 cpl_matrix_delete(coeffs);
01401 cpl_matrix_delete(rhs);
01402 cpl_matrix_delete(result);
01403
01404
01405 for(uj = 0; uj < irplib_hist_get_nbins(hist) - 1; uj++) {
01406 if(irplib_hist_get_value(hist, uj) <= max_interp / 2 &&
01407 irplib_hist_get_value(hist, uj + 1) > max_interp / 2 &&
01408 uj <= maxwhere) {
01409 x1a = uj;
01410 points++;
01411 }
01412 if(irplib_hist_get_value(hist, uj) >= max_interp / 2 &&
01413 irplib_hist_get_value(hist, uj + 1) < max_interp / 2 &&
01414 uj >= maxwhere) {
01415 x2a = uj;
01416 points++;
01417 }
01418 }
01419
01420 if(points != 2) {
01421 cpl_msg_error(cpl_func, "Couldn't compute FWHM");
01422 irplib_hist_delete(hist);
01423 return CPL_ERROR_DATA_NOT_FOUND;
01424 }
01425
01426 x1 = (max_interp / 2 - irplib_hist_get_value(hist, x1a)) /
01427 (irplib_hist_get_value(hist, x1a + 1) -
01428 irplib_hist_get_value(hist, x1a)) + x1a;
01429 x2 = (max_interp / 2 - irplib_hist_get_value(hist, x2a)) /
01430 (irplib_hist_get_value(hist, x2a + 1) -
01431 irplib_hist_get_value(hist, x2a)) + x2a;
01432
01433 *fwhm = (x2 - x1) * irplib_hist_get_bin_size(hist);
01434
01435 *max = max_interp;
01436
01437 *bias = maxwhere_interp * irplib_hist_get_bin_size(hist) +
01438 irplib_hist_get_start(hist);
01439
01440 irplib_hist_delete(hist);
01441
01442 return cpl_error_get_code();
01443 }
01444
01445
01446
01447
01448
01449
01450
01451
01452
01453
01454
01455 static cpl_error_code
01456 irplib_detmon_ronbias_preoverscan(const cpl_imagelist * rawbiases,
01457
01458 cpl_propertylist * qclist,
01459 cpl_image ** synthetic)
01460 {
01461 int i;
01462 int nx, ny;
01463 int nraws;
01464
01465 cpl_vector *meanspre;
01466 cpl_vector *medspre;
01467 cpl_vector *rmsspre;
01468 cpl_vector *meansover;
01469 cpl_vector *medsover;
01470 cpl_vector *rmssover;
01471
01472 cpl_error_code error;
01473
01474 nraws = cpl_imagelist_get_size(rawbiases);
01475 cpl_ensure_code(nraws != -1, CPL_ERROR_ILLEGAL_INPUT);
01476
01477 meanspre = cpl_vector_new(nraws);
01478 medspre = cpl_vector_new(nraws);
01479 rmsspre = cpl_vector_new(nraws);
01480 meansover = cpl_vector_new(nraws);
01481 medsover = cpl_vector_new(nraws);
01482 rmssover = cpl_vector_new(nraws);
01483
01484 nx = cpl_image_get_size_x(cpl_imagelist_get_const(rawbiases, 0));
01485 ny = cpl_image_get_size_y(cpl_imagelist_get_const(rawbiases, 0));
01486 cpl_ensure_code(nx != -1 && ny != -1, CPL_ERROR_ILLEGAL_INPUT);
01487
01488 if(nx < detmon_ronbias_config.prescan_urx ||
01489 nx < detmon_ronbias_config.overscan_urx ||
01490 ny < detmon_ronbias_config.prescan_ury ||
01491 ny < detmon_ronbias_config.overscan_ury) {
01492 cpl_msg_warning(cpl_func, "PREOVERSCAN method not applied. Given "
01493 "limits of prescan and overscan area "
01494 "exceed image size. Please check and rerun.");
01495 return CPL_ERROR_NONE;
01496 }
01497
01498 for(i = 0; i < nraws; i++) {
01499 double mean = 0;
01500 double stdev = 0;
01501
01502 cpl_image *prescan = NULL;
01503 cpl_image *overscan = NULL;
01504
01505 const cpl_image *c_raw = cpl_imagelist_get_const(rawbiases, i);
01506
01507 cpl_ensure_code(c_raw != NULL, CPL_ERROR_ILLEGAL_INPUT);
01508
01509 prescan =
01510 cpl_image_extract(c_raw,
01511 detmon_ronbias_config.prescan_llx,
01512 detmon_ronbias_config.prescan_lly,
01513 detmon_ronbias_config.prescan_urx,
01514 detmon_ronbias_config.prescan_ury);
01515 cpl_ensure_code(prescan != NULL, CPL_ERROR_ILLEGAL_INPUT);
01516 overscan =
01517 cpl_image_extract(c_raw,
01518 detmon_ronbias_config.overscan_llx,
01519 detmon_ronbias_config.overscan_lly,
01520 detmon_ronbias_config.overscan_urx,
01521 detmon_ronbias_config.overscan_ury);
01522 cpl_ensure_code(overscan != NULL, CPL_ERROR_ILLEGAL_INPUT);
01523
01524 if(i == 0) {
01525 *synthetic = irplib_detmon_build_synthetic(prescan, overscan);
01526 cpl_msg_info(cpl_func, "Creating SYNTHETIC frame");
01527 if(*synthetic == NULL) {
01528 cpl_msg_error(cpl_func, "Error creating SYNTHETIC frame");
01529 return CPL_ERROR_UNSPECIFIED;
01530 }
01531 }
01532
01533 error = irplib_get_clean_mean_window(c_raw,
01534 detmon_ronbias_config.
01535 prescan_llx,
01536 detmon_ronbias_config.
01537 prescan_lly,
01538 detmon_ronbias_config.
01539 prescan_urx,
01540 detmon_ronbias_config.
01541 prescan_ury,
01542 detmon_ronbias_config.
01543 stacking_ks_low,
01544 detmon_ronbias_config.
01545 stacking_ks_iter, &mean, &stdev);
01546 cpl_ensure_code(!error, error);
01547
01548 cpl_ensure_code(mean != 0 && stdev != 0, CPL_ERROR_UNSPECIFIED);
01549
01550 error = cpl_vector_set(medspre, i, cpl_image_get_median(prescan));
01551 cpl_ensure_code(!error, error);
01552
01553
01554 error = cpl_vector_set(meanspre, i, mean);
01555 cpl_ensure_code(!error, error);
01556 error = cpl_vector_set(rmsspre, i, stdev);
01557 cpl_ensure_code(!error, error);
01558 error = irplib_get_clean_mean_window(c_raw,
01559 detmon_ronbias_config.
01560 overscan_llx,
01561 detmon_ronbias_config.
01562 overscan_lly,
01563 detmon_ronbias_config.
01564 overscan_urx,
01565 detmon_ronbias_config.
01566 overscan_ury,
01567 detmon_ronbias_config.
01568 stacking_ks_low,
01569 detmon_ronbias_config.
01570 stacking_ks_iter, &mean, &stdev);
01571 cpl_ensure_code(!error, error);
01572
01573 cpl_ensure_code(mean != 0 && stdev != 0, CPL_ERROR_UNSPECIFIED);
01574
01575 error = cpl_vector_set(medsover, i, cpl_image_get_median(overscan));
01576 cpl_ensure_code(!error, error);
01577
01578
01579
01580 error = cpl_vector_set(meansover, i, mean);
01581 cpl_ensure_code(!error, error);
01582 error = cpl_vector_set(rmssover, i, stdev);
01583 cpl_ensure_code(!error, error);
01584
01585 cpl_image_delete(prescan);
01586 cpl_image_delete(overscan);
01587 }
01588
01589 error = cpl_propertylist_append_double(qclist, "ESO QC BIAS PRESCAN MEAN",
01590 cpl_vector_get_mean(meanspre));
01591 cpl_ensure_code(!error, error);
01592 error = cpl_propertylist_append_double(qclist, "ESO QC BIAS PRESCAN MED",
01593 cpl_vector_get_mean(medspre));
01594 cpl_ensure_code(!error, error);
01595 error = cpl_propertylist_append_double(qclist, "ESO QC BIAS PRESCAN RMS",
01596 cpl_vector_get_mean(rmsspre));
01597 cpl_ensure_code(!error, error);
01598
01599 error =
01600 cpl_propertylist_append_double(qclist, "ESO QC BIAS OVERSCAN MEAN",
01601 cpl_vector_get_mean(meansover));
01602 cpl_ensure_code(!error, error);
01603 error = cpl_propertylist_append_double(qclist, "ESO QC BIAS OVERSCAN MED",
01604 cpl_vector_get_mean(medsover));
01605 cpl_ensure_code(!error, error);
01606 error = cpl_propertylist_append_double(qclist, "ESO QC BIAS OVERSCAN RMS",
01607 cpl_vector_get_mean(rmssover));
01608
01609 cpl_ensure_code(!error, error);
01610 error =
01611 cpl_propertylist_append_double(qclist,
01612 "ESO QC BIAS PRESCAN MEAN STDEV",
01613 cpl_vector_get_stdev(meanspre));
01614 cpl_ensure_code(!error, error);
01615 error =
01616 cpl_propertylist_append_double(qclist,
01617 "ESO QC BIAS PRESCAN MED STDEV",
01618 cpl_vector_get_stdev(medspre));
01619 cpl_ensure_code(!error, error);
01620 error =
01621 cpl_propertylist_append_double(qclist,
01622 "ESO QC BIAS PRESCAN RMS STDEV",
01623 cpl_vector_get_stdev(rmsspre));
01624 cpl_ensure_code(!error, error);
01625
01626 error =
01627 cpl_propertylist_append_double(qclist,
01628 "ESO QC BIAS OVERSCAN MEAN STDEV",
01629 cpl_vector_get_stdev(meansover));
01630 cpl_ensure_code(!error, error);
01631 error =
01632 cpl_propertylist_append_double(qclist,
01633 "ESO QC BIAS OVERSCAN MED STDEV",
01634 cpl_vector_get_stdev(medsover));
01635 cpl_ensure_code(!error, error);
01636 error =
01637 cpl_propertylist_append_double(qclist,
01638 "ESO QC BIAS OVERSCAN RMS STDEV",
01639 cpl_vector_get_stdev(rmssover));
01640 cpl_ensure_code(!error, error);
01641
01642 cpl_vector_delete(meanspre);
01643 cpl_vector_delete(medspre);
01644 cpl_vector_delete(rmsspre);
01645 cpl_vector_delete(meansover);
01646 cpl_vector_delete(medsover);
01647 cpl_vector_delete(rmssover);
01648
01649 return CPL_ERROR_NONE;
01650 }
01651
01652
01653
01654
01655
01656
01657
01658
01659
01660
01661
01662
01663 static cpl_error_code
01664 irplib_detmon_ronbias_region(const cpl_imagelist * rawbiases,
01665 const cpl_image * masterbias,
01666 cpl_propertylist * qclist)
01667 {
01668
01669 int nraws = cpl_imagelist_get_size(rawbiases);
01670 int i;
01671
01672 int nx =
01673 cpl_image_get_size_x(cpl_imagelist_get_const(rawbiases, 0));
01674 int ny =
01675 cpl_image_get_size_y(cpl_imagelist_get_const(rawbiases, 0));
01676
01677 cpl_vector *rmssreg;
01678 cpl_error_code error;
01679
01680 const cpl_image * c_raw;
01681 double median, mbias, mstdev;
01682
01683 if(!strcmp(detmon_ronbias_config.pmethod, "DIF")) nraws--;
01684
01685 rmssreg = cpl_vector_new(nraws);
01686
01687 if(nx < detmon_ronbias_config.ref_urx ||
01688 ny < detmon_ronbias_config.ref_ury) {
01689 cpl_msg_warning(cpl_func, "REGION method not applied. Given "
01690 "limits of prescan and overscan area "
01691 "exceed image size. Please check and rerun.");
01692 return CPL_ERROR_NONE;
01693 }
01694
01695 for(i = 0; i < nraws; i++) {
01696 double mean = 0;
01697 double stdev = 0;
01698 if(strcmp(detmon_ronbias_config.pmethod, "DIF")) {
01699 c_raw = cpl_imagelist_get_const(rawbiases, i);
01700 } else {
01701 const cpl_image *c1_raw = cpl_imagelist_get_const(rawbiases, i);
01702 const cpl_image * c2_raw = cpl_imagelist_get_const(rawbiases, i+1);
01703 c_raw = cpl_image_subtract_create(c1_raw, c2_raw);
01704 }
01705 error = irplib_get_clean_mean_window(c_raw,
01706 detmon_ronbias_config.ref_llx,
01707 detmon_ronbias_config.ref_lly,
01708 detmon_ronbias_config.ref_urx,
01709 detmon_ronbias_config.ref_ury,
01710 detmon_ronbias_config.
01711 stacking_ks_low,
01712 detmon_ronbias_config.
01713 stacking_ks_iter, &mean, &stdev);
01714 cpl_ensure_code(!error, error);
01715
01716
01717
01718
01719
01720
01721 error = cpl_vector_set(rmssreg, i, stdev);
01722 cpl_ensure_code(!error, error);
01723 if(!strcmp(detmon_ronbias_config.pmethod, "DIF")) cpl_image_delete((cpl_image *)c_raw);
01724 }
01725
01726 median = cpl_image_get_median_window(masterbias,
01727 detmon_ronbias_config.ref_llx,
01728 detmon_ronbias_config.ref_lly,
01729 detmon_ronbias_config.ref_urx,
01730 detmon_ronbias_config.ref_ury);
01731 error = irplib_get_clean_mean_window(masterbias,
01732 detmon_ronbias_config.ref_llx,
01733 detmon_ronbias_config.ref_lly,
01734 detmon_ronbias_config.ref_urx,
01735 detmon_ronbias_config.ref_ury,
01736 detmon_ronbias_config.
01737 stacking_ks_low,
01738 detmon_ronbias_config.
01739 stacking_ks_iter, &mbias, &mstdev);
01740
01741 error = cpl_propertylist_append_double(qclist, "ESO QC BIAS REGION MED",
01742 median);
01743 cpl_ensure_code(!error, error);
01744
01745 error = cpl_propertylist_append_double(qclist, "ESO QC BIAS REGION BIAS",
01746 mbias);
01747 cpl_ensure_code(!error, error);
01748 error = cpl_propertylist_append_double(qclist, "ESO QC BIAS REGION RMS",
01749 cpl_vector_get_mean(rmssreg));
01750 cpl_ensure_code(!error, error);
01751 error =
01752 cpl_propertylist_append_double(qclist, "ESO QC BIAS REGION RMS STDEV",
01753 cpl_vector_get_stdev(rmssreg));
01754 cpl_ensure_code(!error, error);
01755
01756 cpl_vector_delete(rmssreg);
01757
01758 return cpl_error_get_code();
01759 }
01760
01761
01762
01763
01764
01765
01766
01767
01768
01769
01770
01771
01772 static cpl_image *
01773 irplib_detmon_ronbias_master(const cpl_imagelist * rawbiases,
01774 cpl_propertylist * qclist)
01775 {
01776 double mean = 0;
01777 double stdev = 0;
01778 cpl_image *masterbias = NULL;
01779 cpl_error_code error;
01780
01781 if(!strcmp(detmon_ronbias_config.stacking_method, "MEAN"))
01782 masterbias = cpl_imagelist_collapse_create(rawbiases);
01783 if(!strcmp(detmon_ronbias_config.stacking_method, "MINMAX"))
01784 masterbias =
01785 cpl_imagelist_collapse_minmax_create(rawbiases, 0, 10000);
01786 if(!strcmp(detmon_ronbias_config.stacking_method, "KSIGMA"))
01787 masterbias = cpl_imagelist_collapse_sigclip_create(rawbiases, 3);
01788 if(!strcmp(detmon_ronbias_config.stacking_method, "MEDIAN"))
01789 masterbias = cpl_imagelist_collapse_median_create(rawbiases);
01790
01791 cpl_ensure(masterbias != NULL, CPL_ERROR_DATA_NOT_FOUND, NULL);
01792
01793 error = irplib_get_clean_mean_window(masterbias,
01794 1,
01795 1,
01796 cpl_image_get_size_x(masterbias),
01797 cpl_image_get_size_y(masterbias),
01798 detmon_ronbias_config.
01799 stacking_ks_low,
01800 detmon_ronbias_config.
01801 stacking_ks_iter, &mean, &stdev);
01802 cpl_ensure(!error, error, NULL);
01803
01804 error =
01805 cpl_propertylist_append_double(qclist, "ESO QC MASTER MEAN", mean);
01806
01807 cpl_ensure(!error, error, NULL);
01808 error =
01809 cpl_propertylist_append_double(qclist, "ESO QC MASTER RMS", stdev);
01810
01811 cpl_ensure(!error, error, NULL);
01812
01813 return masterbias;
01814 }
01815
01816
01817
01818
01819
01820
01821
01822
01823
01824
01825
01826
01827 static cpl_error_code
01828 irplib_detmon_ronbias_save(const cpl_parameterlist * parlist,
01829 cpl_frameset * frameset,
01830 const char *recipe_name,
01831 const char *pipeline_name,
01832 const char *procatg_master,
01833 const char *procatg_xstr,
01834 const char *procatg_ystr,
01835 const char *procatg_synth,
01836 const char *package,
01837 const cpl_image * masterbias,
01838 const cpl_image * synthetic,
01839 cpl_propertylist * qclist,
01840 const int flag_sets,
01841 const int which_set, cpl_frameset * usedframes)
01842 {
01843
01844 cpl_frame *ref_frame;
01845 cpl_propertylist *plist;
01846 char *name_o = NULL;
01847
01848
01849 cpl_propertylist *paflist;
01850 cpl_error_code error;
01851
01852
01853
01854
01855
01856
01857 if(!flag_sets) {
01858 name_o = cpl_sprintf("%s_masterbias.fits", recipe_name);
01859 assert(name_o != NULL);
01860 } else {
01861 name_o =
01862 cpl_sprintf("%s_masterbias_set%02d.fits", recipe_name,
01863 which_set);
01864 assert(name_o != NULL);
01865 }
01866
01867
01868 if(cpl_dfs_save_image(frameset, parlist, usedframes, masterbias,
01869 CPL_BPP_IEEE_FLOAT, recipe_name, procatg_master,
01870 qclist, NULL, package, name_o)) {
01871 cpl_msg_error(cpl_func, "Cannot save the product: %s", name_o);
01872
01873 cpl_free(name_o);
01874 cpl_ensure_code(0, CPL_ERROR_FILE_NOT_CREATED);
01875 }
01876
01877
01878 cpl_free(name_o);
01879
01880
01881
01882
01883
01884 if(synthetic) {
01885
01886 if(!flag_sets) {
01887 name_o = cpl_sprintf("%s_synthetic.fits", recipe_name);
01888 assert(name_o != NULL);
01889 } else {
01890 name_o =
01891 cpl_sprintf("%s_synthetic_set%02d.fits", recipe_name,
01892 which_set);
01893 assert(name_o != NULL);
01894 }
01895
01896
01897 if(cpl_dfs_save_image(frameset, parlist, usedframes, synthetic,
01898 CPL_BPP_IEEE_DOUBLE, recipe_name, procatg_synth,
01899 qclist, NULL, package, name_o)) {
01900 cpl_msg_error(cpl_func, "Cannot save the product: %s", name_o);
01901 cpl_free(name_o);
01902
01903 cpl_ensure_code(0, cpl_error_get_code());
01904
01905 }
01906
01907
01908 cpl_free(name_o);
01909 }
01910
01911
01912
01913
01914
01915
01916
01917 ref_frame = cpl_frameset_get_first(frameset);
01918 if((plist = cpl_propertylist_load(cpl_frame_get_filename(ref_frame),
01919 0)) == NULL) {
01920 cpl_msg_error(cpl_func, "getting header from reference frame");
01921 cpl_ensure_code(0, cpl_error_get_code());
01922 }
01923
01924
01925 paflist = cpl_propertylist_new();
01926 cpl_propertylist_copy_property_regexp(paflist, plist,
01927 "^(ARCFILE|MJD-OBS|ESO TPL ID|"
01928 "DATE-OBS|ESO DET DIT|ESO DET NDIT|"
01929 "ESO DET NCORRS|"
01930 "ESO DET MODE NAME)$", 0);
01931
01932 error = cpl_propertylist_append(paflist, qclist);
01933 cpl_ensure_code(!error, error);
01934
01935
01936 if(!flag_sets) {
01937 name_o = cpl_sprintf("%s.paf", recipe_name);
01938 assert(name_o != NULL);
01939 } else {
01940 name_o = cpl_sprintf("%s_set%02d.paf", recipe_name, which_set);
01941 assert(name_o != NULL);
01942 }
01943
01944
01945 if(cpl_dfs_save_paf(pipeline_name, recipe_name, paflist, name_o)) {
01946 cpl_msg_error(cpl_func, "Cannot save the product: %s", name_o);
01947 cpl_free(name_o);
01948 cpl_propertylist_delete(paflist);
01949 cpl_propertylist_delete(plist);
01950 cpl_free(name_o);
01951 cpl_ensure_code(0, CPL_ERROR_FILE_NOT_CREATED);
01952 }
01953
01954 cpl_propertylist_delete(plist);
01955 cpl_propertylist_delete(paflist);
01956 cpl_free(name_o);
01957
01958 if(procatg_xstr == NULL)
01959 cpl_msg_info(cpl_func, "X Structure not created");
01960 if(procatg_ystr == NULL)
01961 cpl_msg_info(cpl_func, "Y Structure not created");
01962
01963 return CPL_ERROR_NONE;
01964 }
01965
01966
01967
01968
01969
01970
01971
01972
01973
01974
01975
01976
01977 int
01978 irplib_detmon_ronbias_dfs_set_groups(cpl_frameset * set, const char *tag)
01979 {
01980 cpl_frame *cur_frame;
01981 const char *cur_tag;
01982 int nframes;
01983 int i;
01984
01985
01986 if(set == NULL)
01987 return -1;
01988
01989
01990 nframes = cpl_frameset_get_size(set);
01991
01992
01993 for(i = 0; i < nframes; i++) {
01994 cur_frame = cpl_frameset_get_frame(set, i);
01995 cur_tag = cpl_frame_get_tag(cur_frame);
01996
01997
01998 if(!strcmp(cur_tag, tag))
01999 cpl_frame_set_group(cur_frame, CPL_FRAME_GROUP_RAW);
02000
02001
02002
02003
02004
02005 }
02006 return 0;
02007 }
02008
02009
02010
02011
02012
02013
02014
02015
02016
02017
02018
02019
02020 static cpl_image *
02021 irplib_detmon_build_synthetic(cpl_image * prescan, cpl_image * overscan)
02022 {
02023 int j;
02024
02025 int distance = detmon_ronbias_config.overscan_urx -
02026 detmon_ronbias_config.prescan_llx + 1;
02027
02028 double *mean_x =
02029 (double *) cpl_malloc(sizeof(double) * distance);
02030
02031 double *xvalues =
02032 (double *) cpl_malloc(sizeof(double) * distance);
02033
02034 cpl_vector *x = NULL;
02035 cpl_vector *y = NULL;
02036
02037 cpl_polynomial *poly = NULL;
02038 cpl_polynomial *poly2 = NULL;
02039
02040 double mse;
02041 int pows[2] = { 0, 0 };
02042
02043 cpl_image *synthetic = NULL;
02044
02045 double initial = 0;
02046
02047
02048 for(j = 0; j < distance; j++) {
02049 *(mean_x + j) = 0;
02050 *(xvalues + j) = j;
02051 }
02052
02053 for(j = 0; j < cpl_image_get_size_x(prescan); j++) {
02054 *(mean_x + j) =
02055 cpl_image_get_mean_window(prescan, j + 1, 1, j + 1,
02056 cpl_image_get_size_y(prescan));
02057 }
02058
02059 for(j = 0; j < cpl_image_get_size_x(overscan); j++) {
02060 *(mean_x + distance - cpl_image_get_size_x(overscan) + j) =
02061 cpl_image_get_mean_window(overscan, j + 1, 1, j + 1,
02062 cpl_image_get_size_y(overscan));
02063 }
02064
02065 x = cpl_vector_wrap(distance, xvalues);
02066 y = cpl_vector_wrap(distance, mean_x);
02067
02068 poly =
02069 cpl_polynomial_fit_1d_create(x, y,
02070 detmon_ronbias_config.preoverscan_degree,
02071 &mse);
02072 cpl_vector_unwrap(x);
02073 cpl_vector_unwrap(y);
02074
02075 initial = *mean_x;
02076
02077 cpl_free(xvalues);
02078 cpl_free(mean_x);
02079
02080 poly2 = cpl_polynomial_new(2);
02081
02082 j = 0;
02083 cpl_polynomial_set_coeff(poly2, pows, cpl_polynomial_get_coeff(poly, &j));
02084
02085 pows[0] = 1;
02086 j = 1;
02087 cpl_polynomial_set_coeff(poly2, pows, cpl_polynomial_get_coeff(poly, &j));
02088
02089 cpl_polynomial_delete(poly);
02090
02091 synthetic =
02092 cpl_image_new(distance, cpl_image_get_size_y(prescan),
02093 CPL_TYPE_DOUBLE);
02094
02095 if(cpl_image_fill_polynomial(synthetic, poly2, initial, 1, 1, 1)) {
02096 cpl_msg_error(cpl_func, "Error creating the synthetic frame");
02097 cpl_polynomial_delete(poly2);
02098 return NULL;
02099 }
02100
02101 cpl_polynomial_delete(poly2);
02102
02103 return synthetic;
02104 }
02105
02106
02107
02108
02109
02110
02111
02112
02113
02114
02115
02116
02117 static cpl_error_code
02118 irplib_detmon_ronbias_dutycycl(const cpl_frameset * frameset,
02119 cpl_propertylist * qclist)
02120 {
02121 const cpl_frame *first;
02122 cpl_propertylist *plistfirst;
02123 double tfirst;
02124 int nraws;
02125 const cpl_frame *last;
02126 cpl_propertylist *plistlast;
02127 double tlast;
02128 double dutycycl;
02129 cpl_error_code error;
02130
02131 first = cpl_frameset_get_first_const(frameset);
02132 cpl_ensure_code(first != NULL, cpl_error_get_code());
02133
02134 plistfirst = cpl_propertylist_load(cpl_frame_get_filename(first), 0);
02135 cpl_ensure_code(plistfirst != NULL, cpl_error_get_code());
02136
02137 tfirst = cpl_propertylist_get_double(plistfirst, "MJD-OBS");
02138 cpl_ensure_code(tfirst != 0, cpl_error_get_code());
02139
02140 nraws = cpl_frameset_get_size(frameset);
02141 cpl_ensure_code(nraws != -1, cpl_error_get_code());
02142
02143 last = cpl_frameset_get_frame_const(frameset, nraws - 1);
02144 cpl_ensure_code(last != NULL, cpl_error_get_code());
02145
02146 plistlast = cpl_propertylist_load(cpl_frame_get_filename(last), 0);
02147 cpl_ensure_code(plistlast != NULL, cpl_error_get_code());
02148
02149 tlast = cpl_propertylist_get_double(plistlast, "MJD-OBS");
02150 cpl_ensure_code(tlast != 0, cpl_error_get_code());
02151
02152 dutycycl = (tlast - tfirst) / (nraws - 1);
02153
02154 error =
02155 cpl_propertylist_append_double(qclist, "ESO QC DUTYCYCL", dutycycl);
02156 cpl_ensure_code(!error, error);
02157
02158 cpl_propertylist_delete(plistfirst);
02159 cpl_propertylist_delete(plistlast);
02160
02161 return CPL_ERROR_NONE;
02162 }
02163
02164
02165
02166
02167
02168
02169
02170
02171
02172
02173
02174
02175
02176
02177 #define HORIZONTAL TRUE
02178
02179 cpl_table *
02180 irplib_detmon_pernoise_reduce(cpl_image * image)
02181 {
02182 int nsamples, nffts;
02183 int i, j;
02184
02185 int status;
02186 float * hanning;
02187 float * data;
02188 float * power;
02189 cpl_image * power_im;
02190 cpl_image * output;
02191 cpl_image * pos_spec;
02192 cpl_table * table;
02193 double freq;
02194 cpl_error_code error;
02195 cpl_image * sub_image;
02196 int nffts_old;
02197
02198 if(detmon_pernoise_config.direction == HORIZONTAL) {
02199 error = cpl_image_flip(image, 1);
02200 cpl_ensure(!error, error, NULL);
02201 }
02202
02203 nsamples = cpl_image_get_size_x(image);
02204 nffts = cpl_image_get_size_y(image);
02205
02206
02207
02208
02209
02210
02211
02212
02213
02214
02215 error = irplib_detmon_pernoise_rm_bg(image, nsamples, nffts);
02216 cpl_ensure(!error, error, NULL);
02217
02218
02219
02220
02221
02222
02223
02224
02225
02226
02227
02228
02229
02230
02231
02232
02233
02234
02235
02236
02237
02238 sub_image = cpl_image_extract(image, nsamples/8 + 1, nffts/8+1,
02239 nsamples*7/8, nffts*7/8);
02240 nffts_old = nffts;
02241 nsamples = cpl_image_get_size_x(sub_image);
02242 nffts = cpl_image_get_size_y(sub_image);
02243
02244
02245
02246
02247
02248
02249
02250
02251
02252
02253
02254
02255
02256
02257
02258
02259 hanning = cpl_malloc(sizeof(float) * nsamples);
02260
02261 for(i = 0; i < nsamples; i++) {
02262 *(hanning + i) = 0.5 - 0.5 * cos(2 * IRPLIB_MATH_PI * (float) i / nsamples);
02263 for(j = 0; j < nffts; j++) {
02264 double value =
02265 cpl_image_get(sub_image, i + 1, j + 1, &status);
02266 error = cpl_image_set(sub_image, i + 1, j + 1, (*(hanning + i)) * value);
02267
02268 cpl_ensure(!error, error, NULL);
02269 }
02270 }
02271
02272 data = cpl_image_get_data_float(sub_image);
02273
02274 power = (float *) cpl_calloc(sizeof(float), nsamples * nffts);
02275
02276
02277
02278 for(i = 0; i < nffts; i++) {
02279 #ifdef HAVE_FFTW
02280 fftwf_complex *fourier =
02281 (fftwf_complex *) fftwf_malloc(sizeof(fftwf_complex) * nsamples);
02282
02283 fftwf_plan p = fftwf_plan_dft_r2c_1d(nsamples,
02284 data + i * nsamples,
02285 fourier,
02286 FFTW_ESTIMATE);
02287 fftwf_execute(p);
02288
02289 for(j = 0; j < nsamples; j++) {
02290 *(power + j + i * nsamples) = fourier[j] * conjf(fourier[j]);
02291 }
02292
02293
02294 fftwf_free(fourier);
02295 #endif
02296 }
02297
02298 power_im = cpl_image_wrap_float(nsamples, nffts, power);
02299 output = cpl_image_collapse_create(power_im, 0);
02300 pos_spec = cpl_image_extract(output, 1, 1, nsamples/2, 1);
02301
02302 cpl_image_unwrap(power_im);
02303 cpl_free(power);
02304
02305 cpl_image_delete(output);
02306
02307 table = cpl_table_new(nsamples/2);
02308 cpl_table_new_column(table, "FREQ", CPL_TYPE_DOUBLE);
02309 cpl_table_new_column(table, "POW", CPL_TYPE_DOUBLE);
02310
02311 freq = detmon_pernoise_config.speed*1000/nffts_old;
02312
02313 for(i = 0; i < nsamples/2; i++) {
02314 error = cpl_table_set(table, "FREQ", i, freq/(nsamples/2)*i);
02315 error = cpl_table_set(table, "POW", i, cpl_image_get(pos_spec, i+1, 1, &status));
02316 }
02317
02318 for(i= 0; i < 5; i++) {
02319 error = cpl_table_set(table, "POW", i, 0.0);
02320 }
02321
02322 cpl_ensure(!error, error, NULL);
02323
02324 cpl_image_delete(pos_spec);
02325
02326 cpl_image_delete(sub_image);
02327 return table;
02328 }
02329 #undef HORIZONTAL
02330
02331
02332
02333
02334
02335
02336
02337
02338
02339
02340
02341
02342
02343
02344
02345
02346 cpl_error_code
02347 irplib_detmon_rm_bpixs(cpl_image ** image,
02348 const double kappa, int nffts, int nsamples)
02349 {
02350 int i, j;
02351
02352 float *data = cpl_image_get_data_float(*image);
02353 int k = 0;
02354 for(i = 0; i < nffts; i++) {
02355 for(j = 0; j < nsamples; j++) {
02356 float neighbours = 0;
02357 int nneighs = 0;
02358 float average = 0;
02359
02360
02361
02362
02363
02364
02365 if(i > 0) {
02366 neighbours += *(data + (i - 1) * nsamples + j);
02367 nneighs++;
02368 }
02369 if(i < nffts - 1) {
02370 neighbours += *(data + (i + 1) * nsamples + j);
02371 nneighs++;
02372 }
02373 if(j > 0) {
02374 neighbours += *(data + i * nsamples + (j - 1));
02375 nneighs++;
02376 }
02377 if(j < nsamples - 1) {
02378 neighbours += *(data + i * nsamples + (j + 1));
02379 nneighs++;
02380 }
02381 average = neighbours / nneighs;
02382 if(average > 0) {
02383 if(*(data + i * nsamples + j) < average * (-1 * kappa) ||
02384 *(data + i * nsamples + j) > average * (kappa)) {
02385 k++;
02386 *(data + i * nsamples + j) = average;
02387 }
02388 }
02389 if(average < 0) {
02390 if(*(data + i * nsamples + j) > average * (-1 * kappa) ||
02391 *(data + i * nsamples + j) < average * (kappa)) {
02392 k++;
02393 *(data + i * nsamples + j) = average;
02394 }
02395 }
02396
02397 }
02398 }
02399
02400
02401 return cpl_error_get_code();
02402
02403 }
02404
02405
02406
02407 #define RECT_RON_HS 4
02408 #define RECT_RON_SAMPLES 100
02409
02410
02411
02412
02413
02414
02415
02416
02417
02418
02419
02420
02421 static cpl_error_code
02422 irplib_flux_get_bias_window(const cpl_image * diff,
02423 const int *zone_def,
02424 int ron_hsize,
02425 int ron_nsamp, double *bias, double *error)
02426 {
02427 const int hsize = ron_hsize < 0 ? RECT_RON_HS : ron_hsize;
02428 const int nsamples =
02429 ron_nsamp < 0 ? RECT_RON_SAMPLES : ron_nsamp;
02430 cpl_bivector *sample_reg;
02431 cpl_vector *rms_list;
02432 int rect[4];
02433 int zone[4];
02434 double *px;
02435 double *py;
02436 double *pr;
02437 int i;
02438
02439
02440 cpl_ensure_code(diff && bias, CPL_ERROR_NULL_INPUT);
02441
02442
02443 if(zone_def != NULL) {
02444 rect[0] = zone_def[0] + hsize + 1;
02445 rect[1] = zone_def[1] - hsize - 1;
02446 rect[2] = zone_def[2] + hsize + 1;
02447 rect[3] = zone_def[3] - hsize - 1;
02448 } else {
02449 rect[0] = hsize + 1;
02450 rect[1] = cpl_image_get_size_x(diff) - hsize - 1;
02451 rect[2] = hsize + 1;
02452 rect[3] = cpl_image_get_size_y(diff) - hsize - 1;
02453 }
02454
02455 cpl_ensure_code(rect[0] < rect[1] && rect[2] < rect[3],
02456 CPL_ERROR_ILLEGAL_INPUT);
02457
02458
02459
02460 sample_reg =
02461 irplib_bivector_gen_rect_poisson(rect, nsamples + 1, nsamples + 1);
02462 cpl_ensure(sample_reg != NULL, CPL_ERROR_ILLEGAL_INPUT,
02463 CPL_ERROR_ILLEGAL_INPUT);
02464
02465 px = cpl_bivector_get_x_data(sample_reg);
02466 py = cpl_bivector_get_y_data(sample_reg);
02467
02468
02469
02470 rms_list = cpl_vector_new(nsamples);
02471 cpl_ensure(rms_list != NULL, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
02472 pr = cpl_vector_get_data(rms_list);
02473
02474 for(i = 0; i < nsamples; i++) {
02475 zone[0] = (int) px[i + 1] - hsize;
02476 zone[1] = (int) px[i + 1] + hsize;
02477 zone[2] = (int) py[i + 1] - hsize;
02478 zone[3] = (int) py[i + 1] + hsize;
02479 pr[i] = cpl_image_get_mean_window(diff,
02480 zone[0], zone[2], zone[1], zone[3]);
02481 }
02482 cpl_bivector_delete(sample_reg);
02483
02484
02485 if(error != NULL)
02486 *error = cpl_vector_get_stdev(rms_list);
02487
02488
02489
02490 *bias = cpl_vector_get_median(rms_list);
02491
02492 cpl_vector_delete(rms_list);
02493
02494 return CPL_ERROR_NONE;
02495 }
02496
02497 #undef RECT_RON_HS
02498 #undef RECT_RON_SAMPLES
02499
02500
02501
02502
02503
02504
02505
02506
02507
02508
02509
02510
02511 static cpl_bivector *
02512 irplib_bivector_gen_rect_poisson(const int *r, const int np, const int homog)
02513 {
02514 double min_dist;
02515 int i;
02516 int gnp;
02517 cpl_bivector *list;
02518 double cand_x, cand_y;
02519 int ok;
02520 int start_ndx;
02521 int xmin, xmax, ymin, ymax;
02522
02523
02524 const int homogc = 0 < homog && homog < np ? homog : np;
02525 double *px;
02526 double *py;
02527
02528
02529 cpl_ensure(r, CPL_ERROR_NULL_INPUT, NULL);
02530 cpl_ensure(np > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
02531
02532 list = cpl_bivector_new(np);
02533 cpl_ensure(list, CPL_ERROR_NULL_INPUT, NULL);
02534 px = cpl_bivector_get_x_data(list);
02535 py = cpl_bivector_get_y_data(list);
02536
02537 xmin = r[0];
02538 xmax = r[1];
02539 ymin = r[2];
02540 ymax = r[3];
02541
02542 min_dist =
02543 IRPLIB_MATH_SQRT1_2 * ((xmax - xmin) * (ymax - ymin) / (double) (homogc + 1));
02544 gnp = 1;
02545 px[0] = 0;
02546 py[0] = 0;
02547
02548
02549 while(gnp < homogc) {
02550
02551 cand_x = cpl_drand() * (xmax - xmin) + xmin;
02552 cand_y = cpl_drand() * (ymax - ymin) + ymin;
02553
02554
02555 ok = 1;
02556 for(i = 0; i < gnp; i++) {
02557 if(pdist(cand_x, cand_y, px[i], py[i]) < min_dist) {
02558
02559 ok = 0;
02560 break;
02561 }
02562 }
02563 if(ok) {
02564
02565 px[gnp] = cand_x;
02566 py[gnp] = cand_y;
02567 gnp++;
02568 }
02569 }
02570
02571
02572
02573 start_ndx = 0;
02574 while(gnp < np) {
02575
02576 cand_x = cpl_drand() * (xmax - xmin) + xmin;
02577 cand_y = cpl_drand() * (ymax - ymin) + ymin;
02578
02579
02580 ok = 1;
02581 for(i = 0; i < homogc; i++) {
02582 if(pdist(cand_x,
02583 cand_y,
02584 px[start_ndx + i], py[start_ndx + i]) < min_dist) {
02585
02586 ok = 0;
02587 break;
02588 }
02589 }
02590 if(ok) {
02591
02592 px[gnp] = cand_x;
02593 py[gnp] = cand_y;
02594 gnp++;
02595 }
02596 }
02597
02598
02599
02600 start_ndx = 0;
02601 while(gnp < np) {
02602
02603 cand_x = cpl_drand() * (xmax - xmin) + xmin;
02604 cand_y = cpl_drand() * (ymax - ymin) + ymin;
02605
02606
02607 ok = 1;
02608 for(i = 0; i < homogc; i++) {
02609 if(pdist(cand_x,
02610 cand_y,
02611 px[start_ndx + i], py[start_ndx + i]) < min_dist) {
02612
02613 ok = 0;
02614 break;
02615 }
02616 }
02617 if(ok) {
02618
02619 px[gnp] = cand_x;
02620 py[gnp] = cand_y;
02621 gnp++;
02622 start_ndx++;
02623 }
02624 }
02625 return list;
02626 }
02627
02628
02629
02630
02631
02632
02633
02634
02635
02636
02637
02638
02639
02640
02641 cpl_error_code
02642 irplib_detmon_pernoise(cpl_frameset * frameset,
02643 const cpl_parameterlist * parlist,
02644 const char * tag,
02645 const char * recipe_name,
02646 const char * pipeline_name,
02647 const char * procatg_tbl,
02648 const char * package,
02649 int (*compare)(const cpl_frame *,
02650 const cpl_frame *))
02651 {
02652 int nsets;
02653 int *selection = NULL;
02654 int i;
02655 cpl_error_code error;
02656
02657 if(irplib_detmon_pernoise_dfs_set_groups(frameset, tag)) {
02658 cpl_msg_error(cpl_func, "Cannot identify RAW and CALIB frames");
02659 }
02660
02661
02662
02663
02664
02665
02666 error = irplib_detmon_pernoise_retrieve_parlist(pipeline_name,
02667 recipe_name, parlist);
02668 cpl_ensure_code(!error, error);
02669
02670
02671 if(compare == NULL)
02672 nsets = 1;
02673 else {
02674 cpl_msg_info(cpl_func, "Identify the different settings");
02675 selection = cpl_frameset_labelise(frameset, compare, &nsets);
02676 if(selection == NULL)
02677 cpl_msg_error(cpl_func, "Cannot labelise input frames");
02678 }
02679
02680 detmon_pernoise_config.nb_extensions = 1;
02681 if(detmon_pernoise_config.exts < 0) {
02682 const cpl_frame *cur_frame =
02683 cpl_frameset_get_first_const(frameset);
02684
02685 detmon_pernoise_config.nb_extensions =
02686 cpl_frame_get_nextensions(cur_frame);
02687 }
02688
02689
02690 for(i = 0; i < nsets; i++) {
02691 int j;
02692 cpl_table ** freq_table;
02693 cpl_propertylist ** qclist =
02694 (cpl_propertylist **)
02695 cpl_malloc(detmon_pernoise_config.nb_extensions *
02696 sizeof(cpl_propertylist *));
02697
02698 cpl_imagelist ** raws = (cpl_imagelist **) cpl_malloc(detmon_pernoise_config.nb_extensions * sizeof(cpl_imagelist *));
02699 cpl_image ** input = (cpl_image **) cpl_malloc(detmon_pernoise_config.nb_extensions * sizeof(cpl_image *));
02700
02701
02702 if(detmon_pernoise_config.mode == 1) {
02703 freq_table =
02704 (cpl_table **) cpl_malloc(detmon_pernoise_config.nb_extensions *
02705 4 * sizeof(cpl_table *));
02706 } else {
02707 freq_table =
02708 (cpl_table **) cpl_malloc(detmon_pernoise_config.nb_extensions *
02709 sizeof(cpl_table *));
02710 }
02711
02712 if(detmon_pernoise_config.exts >= 0) {
02713 *raws =
02714 cpl_imagelist_load_frameset(frameset, CPL_TYPE_FLOAT, 1,
02715 detmon_pernoise_config.exts);
02716 *input = cpl_image_subtract_create(cpl_imagelist_get(*raws,0),
02717 cpl_imagelist_get(*raws,1));
02718 } else {
02719 cpl_imagelist *raws_all_exts =
02720 cpl_imagelist_load_frameset(frameset, CPL_TYPE_FLOAT, 1,
02721 -1);
02722 for(j = 0; j < detmon_pernoise_config.nb_extensions; j++) {
02723 int nframes = cpl_frameset_get_size(frameset);
02724 int k;
02725 for(k = 0; k < nframes; k++) {
02726 cpl_image *image =
02727 cpl_imagelist_unset(raws_all_exts,
02728 (detmon_pernoise_config.
02729 nb_extensions - 1 - j) * k);
02730 cpl_imagelist_set(raws[j], image, k);
02731 }
02732 input[j] =
02733 cpl_image_subtract_create(cpl_imagelist_get(raws[j],0),
02734 cpl_imagelist_get(raws[j],1));
02735 }
02736 }
02737
02738 for(j = 0; j < detmon_pernoise_config.nb_extensions; j++) {
02739 cpl_msg_info(cpl_func, "Starting reduction");
02740 if(detmon_pernoise_config.mode == 1){
02741 int nx = cpl_image_get_size_x(input[j]);
02742 int ny = cpl_image_get_size_y(input[j]);
02743
02744 cpl_image * quad1 = cpl_image_extract(input[j],
02745 1, 1, nx/2, ny/2);
02746 cpl_image * quad2 = cpl_image_extract(input[j],
02747 1, ny/2+1, nx/2, ny);
02748 cpl_image * quad3 = cpl_image_extract(input[j],
02749 nx/2+1, 1, nx, ny/2);
02750 cpl_image * quad4 = cpl_image_extract(input[j],
02751 nx/2+1, ny/2+1, nx, ny);
02752 freq_table[j * 4 ] = irplib_detmon_pernoise_reduce(quad1);
02753 cpl_ensure_code(freq_table[j * 4 ]!= NULL,
02754 cpl_error_get_code());
02755 freq_table[j * 4 + 1] = irplib_detmon_pernoise_reduce(quad2);
02756 cpl_ensure_code(freq_table[j * 4 + 1]!= NULL,
02757 cpl_error_get_code());
02758 freq_table[j * 4 + 2] = irplib_detmon_pernoise_reduce(quad3);
02759 cpl_ensure_code(freq_table[j * 4 + 2]!= NULL,
02760 cpl_error_get_code());
02761 freq_table[j * 4 + 3] = irplib_detmon_pernoise_reduce(quad4);
02762 cpl_ensure_code(freq_table[j * 4 + 3]!= NULL,
02763 cpl_error_get_code());
02764 qclist[j] = cpl_propertylist_new();
02765 error = irplib_detmon_pernoise_qc(qclist[j], freq_table[j], 1);
02766 cpl_ensure_code(!error, error);
02767 error = irplib_detmon_pernoise_qc(qclist[j], freq_table[j+1], 2);
02768 cpl_ensure_code(!error, error);
02769 error = irplib_detmon_pernoise_qc(qclist[j], freq_table[j+2], 3);
02770 cpl_ensure_code(!error, error);
02771 error = irplib_detmon_pernoise_qc(qclist[j], freq_table[j+3], 4);
02772 cpl_ensure_code(!error, error);
02773 } else {
02774 freq_table[j] = irplib_detmon_pernoise_reduce(input[j]);
02775 cpl_ensure_code(freq_table[j]!= NULL, cpl_error_get_code());
02776 qclist[j] = cpl_propertylist_new();
02777 error = irplib_detmon_pernoise_qc(qclist[j], freq_table[j], 0);
02778 cpl_ensure_code(!error, error);
02779 }
02780 }
02781 error = irplib_detmon_pernoise_save(parlist, frameset, recipe_name,
02782 pipeline_name, procatg_tbl,
02783 package, freq_table, qclist, 0,
02784 0, frameset);
02785 cpl_ensure_code(!error, error);
02786
02787 for(j = 0; j < detmon_pernoise_config.nb_extensions; j++) {
02788 cpl_propertylist_delete(qclist[j]);
02789 }
02790 if(detmon_pernoise_config.mode == 1) {
02791 for(j= 0; j < detmon_pernoise_config.nb_extensions * 4; j++) {
02792 cpl_table_delete(freq_table[j]);
02793 }
02794 } else {
02795 for(j= 0; j < detmon_pernoise_config.nb_extensions; j++) {
02796 cpl_table_delete(freq_table[j]);
02797 }
02798 }
02799 }
02800
02801 return cpl_error_get_code();
02802 }
02803
02804
02805
02806
02807
02808
02809
02810
02811
02812
02813
02814
02815 int
02816 irplib_detmon_pernoise_dfs_set_groups(cpl_frameset * set, const char *tag)
02817 {
02818 cpl_frame *cur_frame;
02819 const char *cur_tag;
02820 int nframes;
02821 int i;
02822
02823
02824 if(set == NULL)
02825 return -1;
02826
02827
02828 nframes = cpl_frameset_get_size(set);
02829
02830
02831 for(i = 0; i < nframes; i++) {
02832 cur_frame = cpl_frameset_get_frame(set, i);
02833 cur_tag = cpl_frame_get_tag(cur_frame);
02834
02835
02836 if(!strcmp(cur_tag, tag))
02837 cpl_frame_set_group(cur_frame, CPL_FRAME_GROUP_RAW);
02838
02839
02840
02841
02842
02843 }
02844 return 0;
02845 }
02846
02847
02848
02849
02850
02851
02852
02853
02854
02855
02856
02857
02858 cpl_error_code
02859 irplib_detmon_fill_pernoise_params(cpl_parameterlist * parlist,
02860 const char *recipe_name,
02861 const char *pipeline_name,
02862 int mode,
02863 const char * direction,
02864 double speed,
02865 int llx,
02866 int lly,
02867 int urx,
02868 int ury,
02869 int kappa,
02870 int exts)
02871 {
02872 irplib_detmon_fill_parlist(parlist, recipe_name, pipeline_name, 9,
02873
02874 "mode",
02875 "Mode",
02876 "CPL_TYPE_INT", mode,
02877
02878 "direction",
02879 "Readout direction",
02880 "CPL_TYPE_BOOL", direction,
02881
02882 "speed",
02883 "Readout speed",
02884 "CPL_TYPE_DOUBLE", speed,
02885
02886 "llx",
02887 "(yet unsupported) x coordinate of the lower-left "
02888 "point of the region of interest. If not modified, default value will be 1.",
02889 "CPL_TYPE_INT", llx,
02890 "lly",
02891 "(yet unsupported) y coordinate of the lower-left "
02892 "point of the region of interest. If not modified, default value will be 1.",
02893 "CPL_TYPE_INT", lly,
02894 "urx",
02895 "(yet unsupported) x coordinate of the upper-right "
02896 "point of the region of interest. If not modified, default value will be X dimension of the input image.",
02897 "CPL_TYPE_INT", urx,
02898 "ury",
02899 "(yet unsupported) y coordinate of the upper-right "
02900 "point of the region of interest. If not modified, default value will be Y dimension of the input image.",
02901 "CPL_TYPE_INT", ury,
02902
02903 "kappa",
02904 "Kappa used for determining threshold of bad (hot, cold) pixels",
02905 "CPL_TYPE_INT", kappa,
02906
02907 "exts",
02908 "Activate the multi-exts option",
02909 "CPL_TYPE_INT", exts);
02910
02911 return 0;
02912 }
02913
02914
02915
02916
02917
02918
02919
02920
02921
02922
02923
02924
02925 int
02926 irplib_detmon_fill_pernoise_params_default(cpl_parameterlist * parlist,
02927 const char *recipe_name,
02928 const char *pipeline_name)
02929 {
02930 irplib_detmon_fill_pernoise_params(parlist, recipe_name, pipeline_name,
02931 1,
02932 "CPL_TRUE",
02933 84.5,
02934 -1,
02935 -1,
02936 -1,
02937 -1,
02938 100,
02939 0);
02940
02941 return 0;
02942
02943 }
02944
02945
02946 static cpl_error_code
02947 irplib_detmon_pernoise_retrieve_parlist(const char *pipeline_name,
02948 const char *recipe_name,
02949 const cpl_parameterlist * parlist)
02950 {
02951 char *par_name;
02952 cpl_parameter *par;
02953
02954
02955 detmon_pernoise_config.mode =
02956 irplib_detmon_retrieve_par("mode", pipeline_name, recipe_name,
02957 parlist);
02958
02959
02960 par_name = cpl_sprintf("%s.%s.direction", pipeline_name, recipe_name);
02961 assert(par_name != NULL);
02962 par = cpl_parameterlist_find((cpl_parameterlist *) parlist, par_name);
02963 detmon_pernoise_config.direction = cpl_parameter_get_bool(par);
02964 cpl_free(par_name);
02965
02966
02967 par_name = cpl_sprintf("%s.%s.speed", pipeline_name, recipe_name);
02968 assert(par_name != NULL);
02969 par = cpl_parameterlist_find((cpl_parameterlist *) parlist, par_name);
02970 detmon_pernoise_config.speed = cpl_parameter_get_double(par);
02971 cpl_free(par_name);
02972
02973
02974 detmon_pernoise_config.llx =
02975 irplib_detmon_retrieve_par("llx", pipeline_name, recipe_name,
02976 parlist);
02977
02978
02979 detmon_pernoise_config.lly =
02980 irplib_detmon_retrieve_par("lly", pipeline_name, recipe_name,
02981 parlist);
02982
02983 detmon_pernoise_config.urx =
02984 irplib_detmon_retrieve_par("urx", pipeline_name, recipe_name,
02985 parlist);
02986
02987 detmon_pernoise_config.ury =
02988 irplib_detmon_retrieve_par("ury", pipeline_name, recipe_name,
02989 parlist);
02990
02991 detmon_pernoise_config.kappa =
02992 irplib_detmon_retrieve_par("kappa", pipeline_name, recipe_name,
02993 parlist);
02994
02995
02996 detmon_pernoise_config.exts =
02997 irplib_detmon_retrieve_par("exts", pipeline_name, recipe_name,
02998 parlist);
02999
03000 if(cpl_error_get_code()) {
03001 cpl_msg_error(cpl_func, "Failed to retrieve the input parameters");
03002 cpl_ensure_code(0, CPL_ERROR_DATA_NOT_FOUND);
03003 }
03004
03005 return cpl_error_get_code();
03006 }
03007
03008
03009
03010
03011
03012
03013
03014
03015
03016
03017
03018
03019 static cpl_error_code
03020 irplib_detmon_pernoise_save(const cpl_parameterlist * parlist,
03021 cpl_frameset * frameset,
03022 const char *recipe_name,
03023 const char *pipeline_name,
03024 const char *procatg_tbl,
03025 const char *package,
03026 cpl_table ** freq_table,
03027 cpl_propertylist ** qclist,
03028 const int flag_sets,
03029 const int which_set,
03030 const cpl_frameset * usedframes)
03031 {
03032
03033 cpl_frame *ref_frame;
03034 cpl_propertylist *plist;
03035 char *name_o = NULL;
03036 int i, j;
03037 cpl_propertylist *paflist;
03038 cpl_error_code error;
03039
03040
03041
03042
03043
03044
03045 if(detmon_pernoise_config.mode != 1) {
03046
03047 if(!flag_sets) {
03048 name_o = cpl_sprintf("%s_freq_table.fits", recipe_name);
03049 assert(name_o != NULL);
03050 } else {
03051 name_o =
03052 cpl_sprintf("%s_freq_table_set%02d.fits", recipe_name,
03053 which_set);
03054 assert(name_o != NULL);
03055 }
03056
03057
03058 if(cpl_dfs_save_table(frameset, parlist, usedframes, freq_table[0],
03059 NULL, recipe_name, procatg_tbl, qclist[0], NULL,
03060 package, name_o)) {
03061 cpl_msg_error(cpl_func, "Cannot save the product: %s", name_o);
03062 cpl_free(name_o);
03063 cpl_ensure_code(0, CPL_ERROR_FILE_NOT_CREATED);
03064 }
03065
03066 if(detmon_pernoise_config.exts < 0) {
03067
03068 for(i = 1; i < detmon_pernoise_config.nb_extensions; i++) {
03069 error =
03070 cpl_table_save(freq_table[i], NULL, qclist[i], name_o,
03071 CPL_IO_EXTEND);
03072 cpl_ensure_code(!error, error);
03073 }
03074 }
03075
03076
03077 cpl_free(name_o);
03078
03079 } else {
03080 for (j = 1; j <= 4; j++) {
03081
03082 if(!flag_sets) {
03083 name_o = cpl_sprintf("%s_freq_table_quad%02d.fits",
03084 recipe_name, j);
03085 assert(name_o != NULL);
03086 } else {
03087 name_o =
03088 cpl_sprintf("%s_freq_table_quad%02d_set%02d.fits",
03089 recipe_name, j, which_set);
03090 assert(name_o != NULL);
03091 }
03092
03093
03094 if(cpl_dfs_save_table(frameset, parlist, usedframes,
03095 freq_table[j - 1],
03096 NULL, recipe_name, procatg_tbl, qclist[0], NULL,
03097 package, name_o)) {
03098 cpl_msg_error(cpl_func, "Cannot save the product: %s", name_o);
03099 cpl_free(name_o);
03100 cpl_ensure_code(0, CPL_ERROR_FILE_NOT_CREATED);
03101 }
03102
03103 if(detmon_pernoise_config.exts < 0) {
03104 for(i = 1; i < detmon_pernoise_config.nb_extensions; i++) {
03105 error = cpl_table_save(freq_table[(j-1) + 4 * i],
03106 NULL, qclist[i], name_o,
03107 CPL_IO_EXTEND);
03108 cpl_ensure_code(!error, error);
03109 }
03110 }
03111
03112
03113 cpl_free(name_o);
03114 }
03115
03116 }
03117
03118
03119
03120
03121
03122 ref_frame = cpl_frameset_get_first(frameset);
03123 if((plist = cpl_propertylist_load(cpl_frame_get_filename(ref_frame),
03124 0)) == NULL) {
03125 cpl_msg_error(cpl_func, "getting header from reference frame");
03126 cpl_ensure_code(0, cpl_error_get_code());
03127 }
03128
03129
03130 paflist = cpl_propertylist_new();
03131 cpl_propertylist_copy_property_regexp(paflist, plist,
03132 "^(ARCFILE|MJD-OBS|ESO TPL ID|"
03133 "DATE-OBS|ESO DET DIT|ESO DET NDIT|"
03134 "ESO DET NCORRS|"
03135 "ESO DET MODE NAME)$", 0);
03136
03137 for(i = 0; i < detmon_pernoise_config.nb_extensions; i++) {
03138 cpl_propertylist * c_paflist = cpl_propertylist_duplicate(paflist);
03139 error = cpl_propertylist_append(c_paflist, qclist[i]);
03140 cpl_ensure_code(!error, error);
03141
03142
03143 if(detmon_pernoise_config.exts >= 0) {
03144 if(!flag_sets) {
03145 name_o = cpl_sprintf("%s.paf", recipe_name);
03146 assert(name_o != NULL);
03147 } else {
03148 name_o = cpl_sprintf("%s_set%02d.paf", recipe_name, which_set);
03149 assert(name_o != NULL);
03150 }
03151 } else {
03152 if(!flag_sets) {
03153 name_o = cpl_sprintf("%s_ext%02d.paf", recipe_name, i+1);
03154 assert(name_o != NULL);
03155 } else {
03156 name_o = cpl_sprintf("%s_set%02d_ext%02d.paf", recipe_name, which_set, i+1);
03157 assert(name_o != NULL);
03158 }
03159 }
03160
03161 if(cpl_dfs_save_paf(pipeline_name, recipe_name, c_paflist, name_o)) {
03162 cpl_msg_error(cpl_func, "Cannot save the product: %s", name_o);
03163 cpl_free(name_o);
03164 cpl_propertylist_delete(paflist);
03165 cpl_propertylist_delete(plist);
03166 cpl_free(name_o);
03167 cpl_ensure_code(0, CPL_ERROR_FILE_NOT_CREATED);
03168 }
03169 cpl_propertylist_delete(c_paflist);
03170 cpl_free(name_o);
03171 }
03172
03173 cpl_propertylist_delete(plist);
03174 cpl_propertylist_delete(paflist);
03175
03176 return cpl_error_get_code();
03177 }
03178
03179 static cpl_error_code
03180 irplib_detmon_pernoise_qc(cpl_propertylist * qclist,
03181 cpl_table * table,
03182 int iquad)
03183 {
03184
03185 cpl_error_code error;
03186 char * propname;
03187
03188 double freqs[3] = {0, 0, 0};
03189 double pows[3] = {0, 0, 0};
03190
03191
03192
03193
03194
03195
03196
03197
03198 int nrows = cpl_table_get_nrow(table);
03199 int i;
03200
03201 double * all_freqs = cpl_table_get_data_double(table, "FREQ");
03202 double * all_pows = cpl_table_get_data_double(table, "POW");
03203
03204 for ( i= 1; i< nrows-1; i++){
03205 if (all_pows[i] > pows[0]) {
03206 if(all_pows[i-1] < all_pows[i] && all_pows[i] > all_pows[i+1]){
03207 pows[2]=pows[1];
03208 pows[1]=pows[0];
03209 pows[0]=all_pows[i];
03210
03211 freqs[2]=freqs[1];
03212 freqs[1]=freqs[0];
03213 freqs[0]=all_freqs[i];
03214 }
03215 } else if (all_pows[i] > pows[1]) {
03216 if(all_pows[i-1] < all_pows[i] && all_pows[i] > all_pows[i+1]){
03217 pows[2]=pows[1];
03218 pows[1]=all_pows[i];
03219
03220 freqs[2]=freqs[1];
03221 freqs[1]=all_freqs[i];
03222 }
03223
03224 } else if(all_pows[i] > pows[2]) {
03225 if(all_pows[i-1] < all_pows[i] && all_pows[i] > all_pows[i+1]){
03226 pows[2]=all_pows[i];
03227
03228 freqs[2]=all_freqs[i];
03229 }
03230
03231 }
03232 }
03233
03234 if (detmon_pernoise_config.mode == 1) {
03235 propname = cpl_sprintf("ESO QC FREQ1 %d", iquad);
03236 assert(propname != NULL);
03237 } else {
03238 propname = cpl_sprintf("ESO QC FREQ1");
03239 }
03240
03241 error = cpl_propertylist_append_double(qclist, propname, freqs[0]);
03242 cpl_ensure_code(!error, error);
03243
03244 cpl_free(propname);
03245
03246 if (detmon_pernoise_config.mode == 1) {
03247 propname = cpl_sprintf("ESO QC FREQ2 %d", iquad);
03248 assert(propname != NULL);
03249 } else {
03250 propname = cpl_sprintf("ESO QC FREQ2");
03251 }
03252
03253 error = cpl_propertylist_append_double(qclist, propname, freqs[1]);
03254 cpl_ensure_code(!error, error);
03255
03256 cpl_free(propname);
03257
03258 if (detmon_pernoise_config.mode == 1) {
03259 propname = cpl_sprintf("ESO QC FREQ3 %d", iquad);
03260 assert(propname != NULL);
03261 } else {
03262 propname = cpl_sprintf("ESO QC FREQ3");
03263 }
03264
03265 error = cpl_propertylist_append_double(qclist, propname, freqs[2]);
03266 cpl_ensure_code(!error, error);
03267
03268 cpl_free(propname);
03269
03270 if (detmon_pernoise_config.mode == 1) {
03271 propname = cpl_sprintf("ESO QC POW1 %d", iquad);
03272 assert(propname != NULL);
03273 } else {
03274 propname = cpl_sprintf("ESO QC POW1");
03275 }
03276
03277 error = cpl_propertylist_append_double(qclist, propname, pows[0]);
03278 cpl_ensure_code(!error, error);
03279
03280 cpl_free(propname);
03281
03282 if (detmon_pernoise_config.mode == 1) {
03283 propname = cpl_sprintf("ESO QC POW2 %d", iquad);
03284 assert(propname != NULL);
03285 } else {
03286 propname = cpl_sprintf("ESO QC POW2");
03287 }
03288
03289 error = cpl_propertylist_append_double(qclist, propname, pows[1]);
03290 cpl_ensure_code(!error, error);
03291
03292 cpl_free(propname);
03293
03294 if (detmon_pernoise_config.mode == 1) {
03295 propname = cpl_sprintf("ESO QC POW3 %d", iquad);
03296 assert(propname != NULL);
03297 } else {
03298 propname = cpl_sprintf("ESO QC POW3");
03299 }
03300
03301 error = cpl_propertylist_append_double(qclist, propname, pows[2]);
03302 cpl_ensure_code(!error, error);
03303
03304
03305 cpl_free(propname);
03306
03307 return cpl_error_get_code();
03308 }
03309
03310
03311
03312
03313
03314
03315
03316
03317
03318
03319
03320
03321 cpl_error_code
03322 irplib_detmon_pernoise_rm_bg(cpl_image * image, int nsamples, int nffts)
03323 {
03324 cpl_vector *x = cpl_vector_new(nsamples * nffts);
03325 cpl_vector *y = cpl_vector_new(nsamples * nffts);
03326 cpl_vector *values = cpl_vector_new(nsamples * nffts);
03327
03328 int rejected;
03329 int i, j;
03330 cpl_bivector *xy_pos;
03331 double mse;
03332 cpl_polynomial * poly_2d;
03333 cpl_image * poly_ima;
03334
03335 for(i = 1; i <= nffts; i++) {
03336 for(j = 1; j <= nsamples; j++) {
03337 cpl_vector_set(x, (i - 1) * nsamples + (j - 1), j);
03338 cpl_vector_set(y, (i - 1) * nsamples + (j - 1), i);
03339 cpl_vector_set(values, (i - 1) * nsamples + (j - 1),
03340 cpl_image_get(image, j, i, &rejected));
03341 cpl_ensure_code(!cpl_error_get_code(), cpl_error_get_code());
03342 }
03343 }
03344
03345 xy_pos = cpl_bivector_wrap_vectors(x, y);
03346 poly_2d = cpl_polynomial_fit_2d_create(xy_pos, values, 3, &mse);
03347 poly_ima = cpl_image_new(nsamples, nffts, CPL_TYPE_FLOAT);
03348
03349 cpl_image_fill_polynomial(poly_ima, poly_2d, 1, 1, 1, 1);
03350
03351 cpl_image_subtract(image, poly_ima);
03352
03353
03354
03355
03356
03357
03358
03359
03360
03361
03362
03363
03364
03365
03366
03367
03368 cpl_image_delete(poly_ima);
03369
03370 return cpl_error_get_code();
03371 }
03372
03373
03374
03375
03376
03377
03378
03379
03380
03381
03382
03383
03384 cpl_error_code
03385 irplib_detmon_dark(cpl_frameset * frameset,
03386 const cpl_parameterlist * parlist,
03387 const char * tag,
03388 const char * recipe_name,
03389 const char * pipeline_name,
03390 const char * procatg_master,
03391 const char * procatg_dsnu,
03392 const char * procatg_tbl,
03393 const char * package,
03394 int (*compare)(const cpl_frame *,
03395 const cpl_frame *))
03396 {
03397 int nsets;
03398 int *selection = NULL;
03399 int i;
03400 cpl_error_code error;
03401
03402 if(irplib_detmon_dark_dfs_set_groups(frameset, tag)) {
03403 cpl_msg_error(cpl_func, "Cannot identify RAW and CALIB frames");
03404 }
03405
03406
03407
03408
03409
03410
03411 error = irplib_detmon_retrieve_dark_params(pipeline_name,
03412 recipe_name, parlist);
03413 cpl_ensure_code(!error, error);
03414
03415
03416 if(compare == NULL)
03417 nsets = 1;
03418 else {
03419 cpl_msg_info(cpl_func, "Identify the different settings");
03420 selection = cpl_frameset_labelise(frameset, compare, &nsets);
03421 if(selection == NULL)
03422 cpl_msg_error(cpl_func, "Cannot labelise input frames");
03423 }
03424
03425 detmon_dark_config.nb_extensions = 1;
03426 if(detmon_dark_config.exts < 0) {
03427 const cpl_frame *cur_frame =
03428 cpl_frameset_get_first_const(frameset);
03429
03430 detmon_dark_config.nb_extensions =
03431 cpl_frame_get_nextensions(cur_frame);
03432 }
03433
03434
03435 for(i = 0; i < nsets; i++) {
03436 int *select_dits = NULL;
03437 cpl_frameset *cur_fset =
03438 nsets == 1 ? cpl_frameset_duplicate(frameset) :
03439 cpl_frameset_extract(frameset, selection, i);
03440
03441 int ndits = 0;
03442 int j, k;
03443 cpl_table ** dsnu_table = NULL;
03444 cpl_imagelist ** dsnu = NULL;
03445
03446 cpl_propertylist ** qclist =
03447 (cpl_propertylist **)
03448 cpl_malloc(detmon_dark_config.nb_extensions *
03449 sizeof(cpl_propertylist *));
03450
03451
03452 cpl_imagelist ** masters =
03453 (cpl_imagelist **)
03454 cpl_malloc(detmon_dark_config.nb_extensions *
03455 sizeof(cpl_imagelist *));
03456
03457
03458 if(detmon_dark_config.opt_nir == OPT) {
03459 dsnu_table =
03460 (cpl_table **) cpl_malloc(detmon_dark_config.nb_extensions *
03461 sizeof(cpl_table *));
03462 dsnu =
03463 (cpl_imagelist **)
03464 cpl_malloc(detmon_dark_config.nb_extensions *
03465 sizeof(cpl_imagelist *));
03466 }
03467
03468 select_dits = cpl_frameset_labelise(cur_fset,
03469 irplib_detmon_compare_dits,
03470 &ndits);
03471 cpl_ensure_code(select_dits != NULL, CPL_ERROR_ILLEGAL_INPUT);
03472
03473 if(detmon_dark_config.exts >= 0) {
03474 *masters = cpl_imagelist_new();
03475 if(detmon_dark_config.opt_nir == OPT) {
03476 *dsnu = cpl_imagelist_new();
03477 *dsnu_table = cpl_table_new(ndits);
03478 }
03479 *qclist = cpl_propertylist_new();
03480 cpl_table_new_column(*dsnu_table, "DIT", CPL_TYPE_DOUBLE);
03481 cpl_table_new_column(*dsnu_table, "STDEV", CPL_TYPE_DOUBLE);
03482 } else {
03483 for ( j = 0; j < detmon_dark_config.nb_extensions; j ++) {
03484 masters[j] = cpl_imagelist_new();
03485 if(detmon_dark_config.opt_nir == OPT) {
03486 dsnu[j] = cpl_imagelist_new();
03487 dsnu_table[j] = cpl_table_new(ndits);
03488 }
03489 qclist[j] = cpl_propertylist_new();
03490 cpl_table_new_column(dsnu_table[j], "DIT", CPL_TYPE_DOUBLE);
03491 cpl_table_new_column(dsnu_table[j], "STDEV", CPL_TYPE_DOUBLE);
03492 }
03493 }
03494
03495 for(j = 0; j < ndits; j++) {
03496 cpl_frameset * cur_fdit = cpl_frameset_extract(cur_fset,
03497 select_dits, j);
03498 cpl_imagelist ** raws =
03499 (cpl_imagelist **)
03500 cpl_malloc(detmon_dark_config.nb_extensions *
03501 sizeof(cpl_imagelist *));
03502
03503 if(detmon_dark_config.exts >= 0) {
03504 cpl_image * collapsed;
03505 *raws =
03506 cpl_imagelist_load_frameset(cur_fdit, CPL_TYPE_FLOAT, 1,
03507 detmon_dark_config.exts);
03508 collapsed = cpl_imagelist_collapse_create(*raws);
03509 cpl_imagelist_set(*masters, collapsed, j);
03510 if(detmon_dark_config.opt_nir == OPT) {
03511 irplib_detmon_dark_dsnu(cur_fdit, *dsnu, *dsnu_table,
03512 collapsed, j);
03513 }
03514 irplib_detmon_dark_qc(*qclist, collapsed);
03515 } else {
03516 cpl_imagelist *raws_all_exts =
03517 cpl_imagelist_load_frameset(cur_fdit, CPL_TYPE_FLOAT, 1,
03518 -1);
03519 for(k = 0; k < detmon_dark_config.nb_extensions; k++) {
03520 int nframes = cpl_frameset_get_size(cur_fdit);
03521 int h;
03522 cpl_image * collapsed;
03523 for(h = 0; h < nframes; h++) {
03524 cpl_image *image =
03525 cpl_imagelist_unset(raws_all_exts,
03526 (detmon_dark_config.
03527 nb_extensions - 1 - k) * h);
03528 cpl_imagelist_set(raws[k], image, h);
03529 }
03530 collapsed = cpl_imagelist_collapse_create(raws[k]);
03531 cpl_imagelist_set(masters[k],collapsed, j);
03532 if(detmon_dark_config.opt_nir == OPT) {
03533 irplib_detmon_dark_dsnu(cur_fdit, dsnu[k],
03534 dsnu_table[j], collapsed, j);
03535 }
03536 irplib_detmon_dark_qc(qclist[k], collapsed);
03537 }
03538 }
03539
03540 cpl_frameset_delete(cur_fdit);
03541 for(k = 0; k < detmon_dark_config.nb_extensions; k++) {
03542 cpl_imagelist_delete(raws[k]);
03543 }
03544 cpl_free(raws);
03545 }
03546
03547 cpl_frameset_delete(cur_fset);
03548
03549 irplib_detmon_dark_save(parlist, frameset, recipe_name, pipeline_name,
03550 procatg_master, procatg_tbl, procatg_dsnu,
03551 package, masters, dsnu_table, dsnu, qclist,
03552 0, 0, frameset);
03553
03554 if(detmon_dark_config.opt_nir == OPT) {
03555 for(j = 0; j < detmon_dark_config.nb_extensions; j++) {
03556 cpl_table_delete(dsnu_table[j]);
03557 cpl_imagelist_delete(dsnu[j]);
03558 }
03559 cpl_free(dsnu_table);
03560 cpl_free(dsnu);
03561 }
03562
03563 for(j = 0; j < detmon_dark_config.nb_extensions; j++) {
03564 cpl_propertylist_delete(qclist[j]);
03565 cpl_imagelist_delete(masters[j]);
03566 }
03567 cpl_free(qclist);
03568 cpl_free(masters);
03569 cpl_free(select_dits);
03570
03571 }
03572 cpl_free(selection);
03573
03574 return cpl_error_get_code();
03575 }
03576
03577
03578
03579
03580
03581
03582
03583
03584
03585
03586
03587
03588
03589 int
03590 irplib_detmon_dark_dfs_set_groups(cpl_frameset * set, const char *tag)
03591 {
03592 cpl_frame *cur_frame;
03593 const char *cur_tag;
03594 int nframes;
03595 int i;
03596
03597
03598 if(set == NULL)
03599 return -1;
03600
03601
03602 nframes = cpl_frameset_get_size(set);
03603
03604
03605 for(i = 0; i < nframes; i++) {
03606 cur_frame = cpl_frameset_get_frame(set, i);
03607 cur_tag = cpl_frame_get_tag(cur_frame);
03608
03609
03610 if(!strcmp(cur_tag, tag))
03611 cpl_frame_set_group(cur_frame, CPL_FRAME_GROUP_RAW);
03612
03613
03614
03615
03616
03617 }
03618 return 0;
03619 }
03620
03621
03622
03623
03624
03625
03626
03627
03628
03629
03630
03631
03632 static cpl_error_code
03633 irplib_detmon_retrieve_dark_params(const char *pipeline_name,
03634 const char *recipe_name,
03635 const cpl_parameterlist * parlist)
03636 {
03637 char *par_name;
03638 cpl_parameter *par;
03639
03640
03641 par_name = cpl_sprintf("%s.%s.ron.method", pipeline_name, recipe_name);
03642 assert(par_name != NULL);
03643 par = cpl_parameterlist_find((cpl_parameterlist *) parlist, par_name);
03644 detmon_dark_config.ron_method = cpl_parameter_get_string(par);
03645 cpl_free(par_name);
03646
03647
03648 par_name = cpl_sprintf("%s.%s.dsnu.method", pipeline_name, recipe_name);
03649 assert(par_name != NULL);
03650 par = cpl_parameterlist_find((cpl_parameterlist *) parlist, par_name);
03651 detmon_dark_config.dsnu_method = cpl_parameter_get_string(par);
03652 cpl_free(par_name);
03653
03654
03655 par_name = cpl_sprintf("%s.%s.opt_nir", pipeline_name, recipe_name);
03656 assert(par_name != NULL);
03657 par = cpl_parameterlist_find((cpl_parameterlist *) parlist, par_name);
03658 detmon_dark_config.opt_nir = cpl_parameter_get_bool(par);
03659 cpl_free(par_name);
03660
03661
03662 detmon_dark_config.exts =
03663 irplib_detmon_retrieve_par("exts", pipeline_name, recipe_name,
03664 parlist);
03665
03666 if(cpl_error_get_code()) {
03667 cpl_msg_error(cpl_func, "Failed to retrieve the input parameters");
03668 cpl_ensure_code(0, CPL_ERROR_DATA_NOT_FOUND);
03669 }
03670
03671
03672 return CPL_ERROR_NONE;
03673 }
03674
03675
03676
03677
03678
03679
03680
03681
03682
03683
03684
03685
03686
03687 cpl_error_code
03688 irplib_detmon_fill_dark_params(cpl_parameterlist * parlist,
03689 const char *recipe_name,
03690 const char *pipeline_name,
03691 const char * ron_method,
03692 const char * dsnu_method,
03693 const char * opt_nir,
03694 int exts)
03695 {
03696 irplib_detmon_fill_parlist(parlist, recipe_name, pipeline_name, 4,
03697
03698 "ron.method",
03699 "RON method",
03700 "CPL_TYPE_STRING", ron_method,
03701
03702 "dsnu.method",
03703 "DSNU method",
03704 "CPL_TYPE_STRING", dsnu_method,
03705
03706 "opt_nir",
03707 "OPT or NIR",
03708 "CPL_TYPE_BOOL", opt_nir,
03709
03710 "exts",
03711 "Activate the multi-exts option",
03712 "CPL_TYPE_INT", exts);
03713
03714 return cpl_error_get_code();
03715 }
03716
03717
03718
03719
03720
03721
03722
03723
03724
03725
03726
03727
03728 int
03729 irplib_detmon_fill_dark_params_default(cpl_parameterlist * parlist,
03730 const char *recipe_name,
03731 const char *pipeline_name)
03732 {
03733 irplib_detmon_fill_dark_params(parlist, recipe_name, pipeline_name,
03734 "SIMPLE",
03735 "STDEV",
03736 "CPL_FALSE",
03737 0);
03738 return cpl_error_get_code();
03739 }
03740
03741
03742
03743
03744
03745
03746
03747
03748
03749
03750
03751
03752 cpl_error_code
03753 irplib_detmon_dark_dsnu(cpl_frameset * cur_fdit,
03754 cpl_imagelist * dsnu,
03755 cpl_table * dsnu_table,
03756 cpl_image * collapsed,
03757 int pos)
03758 {
03759 cpl_frame * first = cpl_frameset_get_first(cur_fdit);
03760 cpl_propertylist * plist =
03761 cpl_propertylist_load(cpl_frame_get_filename(first), 0);
03762 double dit = irplib_pfits_get_exptime(plist);
03763 double mean = cpl_image_get_mean(collapsed);
03764
03765 cpl_image * dsnu_map =
03766 cpl_image_subtract_scalar_create(collapsed, mean);
03767 double stdev;
03768 cpl_image_divide_scalar(dsnu_map, mean);
03769 stdev = cpl_image_get_stdev(dsnu_map);
03770
03771 cpl_imagelist_set(dsnu, dsnu_map, pos);
03772
03773 cpl_table_set(dsnu_table, "DIT", pos, dit);
03774 cpl_table_set(dsnu_table, "STDEV", pos, stdev);
03775
03776 cpl_propertylist_delete(plist);
03777
03778 return cpl_error_get_code();
03779
03780 }
03781
03782
03783
03784
03785
03786
03787
03788
03789
03790
03791
03792
03793 static cpl_error_code
03794 irplib_detmon_dark_save(const cpl_parameterlist * parlist,
03795 cpl_frameset * frameset,
03796 const char *recipe_name,
03797 const char *pipeline_name,
03798 const char *procatg_master,
03799 const char *procatg_tbl,
03800 const char *procatg_dsnu,
03801 const char *package,
03802 cpl_imagelist ** masters,
03803 cpl_table ** dsnu_table,
03804 cpl_imagelist ** dsnu,
03805 cpl_propertylist ** qclist,
03806 const int flag_sets,
03807 const int which_set,
03808 const cpl_frameset * usedframes)
03809 {
03810
03811 cpl_frame *ref_frame;
03812 cpl_propertylist *plist;
03813 char *name_o = NULL;
03814 int i, j;
03815 cpl_propertylist *paflist;
03816 cpl_error_code error;
03817 int nb_images;
03818
03819
03820
03821
03822
03823 nb_images = cpl_imagelist_get_size(masters[0]);
03824 cpl_ensure_code(nb_images > 0, CPL_ERROR_DATA_NOT_FOUND);
03825
03826
03827 for(i = 0; i < nb_images; i++) {
03828
03829 if(!flag_sets) {
03830 name_o =
03831 cpl_sprintf("%s_master_dit_%d.fits", recipe_name, i+1);
03832 assert(name_o != NULL);
03833 } else {
03834 name_o =
03835 cpl_sprintf("%s_master_dit_%d_set%02d.fits",
03836 recipe_name, i, which_set);
03837 assert(name_o != NULL);
03838 }
03839
03840
03841
03842 if(detmon_dark_config.exts >= 0) {
03843 if(cpl_dfs_save_image
03844 (frameset, parlist, usedframes,
03845 cpl_imagelist_get(*masters, i), CPL_BPP_IEEE_FLOAT,
03846 recipe_name, procatg_master, qclist[0], NULL, package,
03847 name_o)) {
03848 cpl_msg_error(cpl_func, "Cannot save the product: %s",
03849 name_o);
03850 cpl_free(name_o);
03851 cpl_ensure_code(0, CPL_ERROR_FILE_NOT_CREATED);
03852
03853 }
03854 } else {
03855 if(cpl_dfs_save_image(frameset, parlist, usedframes, NULL,
03856 CPL_BPP_IEEE_FLOAT, recipe_name,
03857 procatg_master, qclist[0], NULL,
03858 package, name_o)) {
03859 cpl_msg_error(cpl_func, "Cannot save the product: %s",
03860 name_o);
03861 cpl_free(name_o);
03862 cpl_ensure_code(0, CPL_ERROR_FILE_NOT_CREATED);
03863 }
03864
03865 for(j = 0; j < detmon_dark_config.nb_extensions; j++) {
03866 error =
03867 cpl_image_save(cpl_imagelist_get(masters[j], i),
03868 name_o, CPL_BPP_IEEE_FLOAT, qclist[j],
03869 CPL_IO_EXTEND);
03870 cpl_ensure_code(!error, error);
03871 }
03872 }
03873 }
03874
03875 cpl_free(name_o);
03876
03877 if (detmon_dark_config.opt_nir == OPT) {
03878
03879
03880
03881
03882
03883
03884 if(!flag_sets) {
03885 name_o = cpl_sprintf("%s_dsnu_table.fits", recipe_name);
03886 assert(name_o != NULL);
03887 } else {
03888 name_o =
03889 cpl_sprintf("%s_dsnu_table_set%02d.fits", recipe_name,
03890 which_set);
03891 assert(name_o != NULL);
03892 }
03893
03894
03895 if(cpl_dfs_save_table(frameset, parlist, usedframes, dsnu_table[0],
03896 NULL, recipe_name, procatg_tbl, qclist[0], NULL,
03897 package, name_o)) {
03898 cpl_msg_error(cpl_func, "Cannot save the product: %s", name_o);
03899 cpl_free(name_o);
03900 cpl_ensure_code(0, CPL_ERROR_FILE_NOT_CREATED);
03901 }
03902
03903 if(detmon_dark_config.exts < 0) {
03904
03905 for(i = 1; i < detmon_dark_config.nb_extensions; i++) {
03906 error =
03907 cpl_table_save(dsnu_table[i], NULL, qclist[i], name_o,
03908 CPL_IO_EXTEND);
03909 cpl_ensure_code(!error, error);
03910 }
03911 }
03912
03913
03914 cpl_free(name_o);
03915
03916
03917
03918
03919
03920 for(i = 0; i < nb_images; i++) {
03921
03922 if(!flag_sets) {
03923 name_o =
03924 cpl_sprintf("%s_dsnu_map_dit_%d.fits", recipe_name, i+1);
03925 assert(name_o != NULL);
03926 } else {
03927 name_o =
03928 cpl_sprintf("%s_dsnu_map_dit_%d_set%02d.fits",
03929 recipe_name, i, which_set);
03930 assert(name_o != NULL);
03931 }
03932
03933
03934
03935 if(detmon_dark_config.exts >= 0) {
03936 if(cpl_dfs_save_image
03937 (frameset, parlist, usedframes,
03938 cpl_imagelist_get(*dsnu, i), CPL_BPP_IEEE_FLOAT,
03939 recipe_name, procatg_dsnu, qclist[0], NULL, package,
03940 name_o)) {
03941 cpl_msg_error(cpl_func, "Cannot save the product: %s",
03942 name_o);
03943 cpl_free(name_o);
03944 cpl_ensure_code(0, CPL_ERROR_FILE_NOT_CREATED);
03945
03946 }
03947 } else {
03948 if(cpl_dfs_save_image(frameset, parlist, usedframes, NULL,
03949 CPL_BPP_IEEE_FLOAT, recipe_name,
03950 procatg_dsnu, qclist[0], NULL,
03951 package, name_o)) {
03952 cpl_msg_error(cpl_func, "Cannot save the product: %s",
03953 name_o);
03954 cpl_free(name_o);
03955 cpl_ensure_code(0, CPL_ERROR_FILE_NOT_CREATED);
03956 }
03957
03958 for(j = 0; j < detmon_dark_config.nb_extensions; j++) {
03959 error =
03960 cpl_image_save(cpl_imagelist_get(dsnu[j], i),
03961 name_o, CPL_BPP_IEEE_FLOAT, qclist[j],
03962 CPL_IO_EXTEND);
03963 cpl_ensure_code(!error, error);
03964 }
03965 }
03966 }
03967
03968 cpl_free(name_o);
03969
03970 }
03971
03972
03973
03974
03975
03976
03977 ref_frame = cpl_frameset_get_first(frameset);
03978 if((plist = cpl_propertylist_load(cpl_frame_get_filename(ref_frame),
03979 0)) == NULL) {
03980 cpl_msg_error(cpl_func, "getting header from reference frame");
03981 cpl_ensure_code(0, cpl_error_get_code());
03982 }
03983
03984
03985 paflist = cpl_propertylist_new();
03986 cpl_propertylist_copy_property_regexp(paflist, plist,
03987 "^(ARCFILE|MJD-OBS|ESO TPL ID|"
03988 "DATE-OBS|ESO DET DIT|ESO DET NDIT|"
03989 "ESO DET NCORRS|"
03990 "ESO DET MODE NAME)$", 0);
03991
03992 for(i = 0; i < detmon_dark_config.nb_extensions; i++) {
03993 cpl_propertylist * c_paflist = cpl_propertylist_duplicate(paflist);
03994 error = cpl_propertylist_append(c_paflist, qclist[i]);
03995 cpl_ensure_code(!error, error);
03996
03997
03998 if(detmon_dark_config.exts >= 0) {
03999 if(!flag_sets) {
04000 name_o = cpl_sprintf("%s.paf", recipe_name);
04001 assert(name_o != NULL);
04002 } else {
04003 name_o = cpl_sprintf("%s_set%02d.paf", recipe_name, which_set);
04004 assert(name_o != NULL);
04005 }
04006 } else {
04007 if(!flag_sets) {
04008 name_o = cpl_sprintf("%s_ext%02d.paf", recipe_name, i+1);
04009 assert(name_o != NULL);
04010 } else {
04011 name_o = cpl_sprintf("%s_set%02d_ext%02d.paf", recipe_name, which_set, i+1);
04012 assert(name_o != NULL);
04013 }
04014 }
04015
04016 if(cpl_dfs_save_paf(pipeline_name, recipe_name, c_paflist, name_o)) {
04017 cpl_msg_error(cpl_func, "Cannot save the product: %s", name_o);
04018 cpl_free(name_o);
04019 cpl_propertylist_delete(paflist);
04020 cpl_propertylist_delete(plist);
04021 cpl_free(name_o);
04022 cpl_ensure_code(0, CPL_ERROR_FILE_NOT_CREATED);
04023 }
04024 cpl_propertylist_delete(c_paflist);
04025 cpl_free(name_o);
04026 }
04027
04028 cpl_propertylist_delete(plist);
04029 cpl_propertylist_delete(paflist);
04030
04031 return cpl_error_get_code();
04032 }
04033
04034 cpl_error_code
04035 irplib_detmon_dark_qc(cpl_propertylist * qclist,
04036 cpl_image * collapsed)
04037 {
04038 double mean = cpl_image_get_mean(collapsed);
04039 double stdev = cpl_image_get_stdev(collapsed);
04040
04041 cpl_error_code error;
04042
04043 error = cpl_propertylist_append_double(qclist, "ESO QC DARK", mean);
04044 cpl_ensure_code(!error, error);
04045
04046 error = cpl_propertylist_append_double(qclist, "ESO QC STDEV", stdev);
04047 cpl_ensure_code(!error, error);
04048
04049 return cpl_error_get_code();
04050 }
04051