27#if !defined(_XOPEN_SOURCE) || (_XOPEN_SOURCE - 0) < 500
28#define _XOPEN_SOURCE 500
31#include "hdrl_utils.h"
32#include "hdrl_types.h"
33#include "hdrl_elemop.h"
34#include "hdrl_prototyping.h"
35#include "hdrl_imagelist.h"
70 static const char * hdrl_license =
71 "This file is part of the HDRL Instrument Pipeline\n"
72 "Copyright (C) 2012,2013 European Southern Observatory\n"
74 "This program is free software; you can redistribute it and/or modify\n"
75 "it under the terms of the GNU General Public License as published by\n"
76 "the Free Software Foundation; either version 2 of the License, or\n"
77 "(at your option) any later version.\n"
79 "This program is distributed in the hope that it will be useful,\n"
80 "but WITHOUT ANY WARRANTY; without even the implied warranty of\n"
81 "MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the\n"
82 "GNU General Public License for more details.\n"
84 "You should have received a copy of the GNU General Public License\n"
85 "along with this program; if not, write to the Free Software\n"
86 "Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, \n"
92#define EOP_LINE_SIZE 188
114 cpl_ensure (eop_data, CPL_ERROR_NULL_INPUT, NULL);
116 if(!((data_length) % EOP_LINE_SIZE == 0))
118 cpl_error_set_message(cpl_func, CPL_ERROR_NULL_INPUT,
119 "Raw data doesn't have a fixed record width");
124 cpl_size n_entries = (data_length) / EOP_LINE_SIZE;
125 cpl_table * eop_table = cpl_table_new (n_entries);
126 cpl_msg_info (cpl_func,
" EOP data has a total of %"CPL_SIZE_FORMAT
" entries", n_entries);
129 cpl_table_new_column (eop_table,
"MJD", CPL_TYPE_DOUBLE);
130 cpl_table_new_column (eop_table,
"PMX", CPL_TYPE_DOUBLE);
131 cpl_table_new_column (eop_table,
"PMY", CPL_TYPE_DOUBLE);
132 cpl_table_new_column (eop_table,
"DUT", CPL_TYPE_DOUBLE);
133 cpl_table_new_column (eop_table,
"FLAG", CPL_TYPE_STRING);
136 cpl_table_set_column_unit (eop_table,
"MJD",
"d");
137 cpl_table_set_column_unit (eop_table,
"PMX",
"arcsec");
138 cpl_table_set_column_unit (eop_table,
"PMY",
"arcsec");
139 cpl_table_set_column_unit (eop_table,
"DUT",
"s");
142 for(cpl_size i=0; i<n_entries; i++)
145 strncpy(flag, eop_data+i*EOP_LINE_SIZE+16, 1);
147 cpl_table_set_string(eop_table,
"FLAG", i, flag);
149 cpl_table_set_double(eop_table,
"MJD", i, atof(eop_data+i*EOP_LINE_SIZE+7));
150 if(!strncmp(flag,
"I", 1) || !strncmp(flag,
"P", 1))
152 cpl_table_set_double(eop_table,
"PMX", i, atof(eop_data+i*EOP_LINE_SIZE+18));
153 cpl_table_set_double(eop_table,
"PMY", i, atof(eop_data+i*EOP_LINE_SIZE+37));
154 cpl_table_set_double(eop_table,
"DUT", i, atof(eop_data+i*EOP_LINE_SIZE+58));
159 cpl_table_unselect_all (eop_table);
160 cpl_table_or_selected_invalid (eop_table,
"PMX");
161 cpl_table_or_selected_invalid (eop_table,
"PMY");
162 cpl_table_or_selected_invalid (eop_table,
"DUT");
163 cpl_msg_info (cpl_func,
"Found %lld invalid", cpl_table_count_selected (eop_table));
164 cpl_table_erase_selected (eop_table);
181} hdrl_rect_region_parameter;
184static hdrl_parameter_typeobj hdrl_rect_region_parameter_type = {
185 HDRL_PARAMETER_RECT_REGION,
186 (hdrl_alloc *)&cpl_malloc,
187 (hdrl_free *)&cpl_free,
189 sizeof(hdrl_rect_region_parameter),
207hdrl_maglim_kernel_create(
const cpl_size kernel_sx,
const cpl_size kernel_sy,
210 cpl_ensure(kernel_sx > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
211 cpl_ensure(kernel_sy > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
212 cpl_ensure(fwhm > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
214 cpl_matrix * kernel = cpl_matrix_new(kernel_sx, kernel_sy);
216 double* pkernel = cpl_matrix_get_data(kernel);
222 double sigma_to_fwhm = sqrt( 4. * log(4) );
224 double factor = (fwhm / sigma_to_fwhm);
226 factor = 2.0 * factor ;
232 for(cpl_size j = 0; j < kernel_sy; j++) {
233 y = j * dy - 0.5 * (kernel_sy - 1);
235 for(cpl_size i = 0; i < kernel_sx; i++) {
236 x = i * dx - 0.5 * (kernel_sx - 1);
238 arg_exp = (x * x + y * y) / factor;
239 pkernel[kernel_sx * j + i] = exp( -arg_exp );
267hdrl_extend_image(
const cpl_image* image,
const cpl_size border_nx,
268 const cpl_size border_ny,
269 const hdrl_image_extend_method image_extend_method)
271 cpl_ensure(image != NULL, CPL_ERROR_NULL_INPUT, NULL);
272 cpl_ensure(border_nx > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
273 cpl_ensure(border_ny > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
275 cpl_ensure(image_extend_method == HDRL_IMAGE_EXTEND_NEAREST ||
276 image_extend_method == HDRL_IMAGE_EXTEND_MIRROR, CPL_ERROR_ILLEGAL_INPUT, NULL);
278 cpl_type type = cpl_image_get_type(image);
279 cpl_size input_sx = cpl_image_get_size_x(image);
280 cpl_size input_sy = cpl_image_get_size_y(image);
282 cpl_ensure(2 * border_nx <= input_sx, CPL_ERROR_ILLEGAL_INPUT, NULL);
283 cpl_ensure(2 * border_ny <= input_sy, CPL_ERROR_ILLEGAL_INPUT, NULL);
285 cpl_size output_sx = input_sx + 2 * border_nx;
286 cpl_size output_sy = input_sy + 2 * border_ny;
288 cpl_msg_debug(cpl_func,
"Extend image:");
289 cpl_msg_indent_more();
290 cpl_msg_debug(cpl_func,
"Border sizes (x, y): (%lld, %lld)", border_nx, border_ny);
291 cpl_msg_debug(cpl_func,
"Input image (x, y): (%lld, %lld)", input_sx, input_sy);
292 cpl_msg_debug(cpl_func,
"Output image (x, y): (%lld, %lld)", output_sx, output_sy);
293 cpl_msg_indent_less();
295 cpl_image* output_image = cpl_image_new(output_sx, output_sy, type);
299 if ( image_extend_method == HDRL_IMAGE_EXTEND_NEAREST ) {
302 cpl_image_copy(output_image, image, border_nx + 1, border_ny + 1 );
305 for(cpl_size j = 1; j <= input_sy; j++) {
307 value = cpl_image_get(image, 1, j, &status);
308 cpl_image_fill_window(output_image, 1, border_ny +j,
309 border_nx, border_ny + j, value);
311 value = cpl_image_get(image, input_sx, j, &status);
312 cpl_image_fill_window(output_image, input_sx + border_nx, border_ny +j,
313 output_sx, border_ny + j, value);
317 for(cpl_size i = 1; i <= output_sx; i++) {
319 value = cpl_image_get(output_image, i, border_ny + 1, &status);
320 cpl_image_fill_window(output_image, i, 1, i, border_ny, value);
323 value = cpl_image_get(output_image, i, border_ny + input_sy, &status);
324 cpl_image_fill_window(output_image, i, border_ny + input_sy, i,
328 }
else if ( image_extend_method == HDRL_IMAGE_EXTEND_MIRROR ) {
331 cpl_image_copy(output_image, image, border_nx + 1, border_ny + 1);
334 cpl_image* image_ex = cpl_image_extract(image, 1, 1, border_nx, input_sy);
335 cpl_image_flip(image_ex, 2);
336 cpl_image_copy(output_image, image_ex, 1, border_ny + 1);
337 cpl_image_delete(image_ex);
340 image_ex = cpl_image_extract(image, input_sx- border_nx + 1, 1, input_sx,
342 cpl_image_flip(image_ex, 2);
343 cpl_image_copy(output_image, image_ex, input_sx + border_nx + 1,
345 cpl_image_delete(image_ex);
348 image_ex = cpl_image_extract(output_image, 1, input_sy, output_sx,
349 input_sy + border_ny);
350 cpl_image_flip(image_ex, 0);
351 cpl_image_copy(output_image, image_ex, 1, input_sy + border_ny + 1);
352 cpl_image_delete(image_ex);
355 image_ex = cpl_image_extract(output_image, 1, border_ny + 1, output_sx,
357 cpl_image_flip(image_ex, 0);
358 cpl_image_copy(output_image, image_ex, 1, 1);
359 cpl_image_delete(image_ex);
385hdrl_image_convolve(
const cpl_image * input_image,
const cpl_matrix * kernel,
386 const hdrl_image_extend_method image_extend_method)
388 cpl_image * expanded_image;
389 cpl_image * filtered_image;
390 cpl_image * expanded_filtered_image;
393 cpl_ensure(input_image != NULL, CPL_ERROR_NULL_INPUT, NULL);
394 cpl_ensure(kernel != NULL, CPL_ERROR_NULL_INPUT, NULL);
395 cpl_ensure(image_extend_method == HDRL_IMAGE_EXTEND_NEAREST ||
396 image_extend_method == HDRL_IMAGE_EXTEND_MIRROR, CPL_ERROR_ILLEGAL_INPUT, NULL);
398 cpl_size kernel_nx = cpl_matrix_get_ncol(kernel);
399 cpl_size kernel_ny = cpl_matrix_get_nrow(kernel);
400 cpl_ensure(kernel_nx >= 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
401 cpl_ensure(kernel_ny >= 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
404 cpl_ensure((kernel_nx&1) == 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
405 cpl_ensure((kernel_ny&1) == 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
407 const int border_nx = (int)((kernel_nx - 1)/2.);
408 const int border_ny = (int)((kernel_ny - 1)/2.);
411 hdrl_extend_image(input_image, border_nx, border_ny, image_extend_method);
415 expanded_filtered_image =
416 hdrl_parallel_filter_image(expanded_image, kernel, NULL,
420 if(expanded_filtered_image == NULL) {
421 cpl_image_delete(expanded_filtered_image);
422 cpl_image_delete(expanded_image);
427 cpl_msg_debug(cpl_func,
"Extract original image from expanded mask, region "\
428 "[llx: %lld, lly: %lld, urx: %lld, ury: %lld",
429 kernel_nx+1, kernel_ny + 1,
430 cpl_image_get_size_x(input_image) + kernel_nx,
431 cpl_image_get_size_y(input_image) + kernel_ny);
433 filtered_image = cpl_image_extract(expanded_filtered_image, border_nx + 1,
435 cpl_image_get_size_x(input_image) + border_nx,
436 cpl_image_get_size_y(input_image) + border_ny);
439 cpl_image_delete(expanded_filtered_image);
440 cpl_image_delete(expanded_image);
442 return filtered_image;
465 hdrl_rect_region_parameter * p = (hdrl_rect_region_parameter *)
466 hdrl_parameter_new(&hdrl_rect_region_parameter_type);
471 return (hdrl_parameter *)p;
486 hdrl_parameter * rect_region,
492 hdrl_rect_region_parameter * p = (hdrl_rect_region_parameter *)rect_region ;
508 return hdrl_parameter_check_type(self, &hdrl_rect_region_parameter_type);
520 cpl_ensure(p, CPL_ERROR_NULL_INPUT, -1) ;
521 return ((
const hdrl_rect_region_parameter *)p)->llx;
529 cpl_ensure(p, CPL_ERROR_NULL_INPUT, -1) ;
530 return ((
const hdrl_rect_region_parameter *)p)->lly;
538 cpl_ensure(p, CPL_ERROR_NULL_INPUT, -1) ;
539 return ((
const hdrl_rect_region_parameter *)p)->urx;
547 cpl_ensure(p, CPL_ERROR_NULL_INPUT, -1) ;
548 return ((
const hdrl_rect_region_parameter *)p)->ury;
561 const hdrl_parameter * param,
562 const cpl_size max_x,
563 const cpl_size max_y)
565 const hdrl_rect_region_parameter * param_loc =
566 (
const hdrl_rect_region_parameter *)param ;
568 cpl_error_ensure(param != NULL, CPL_ERROR_NULL_INPUT,
569 return CPL_ERROR_NULL_INPUT,
"NULL Input Parameters");
571 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
572 "Expected Rect Region parameter") ;
573 cpl_error_ensure(param_loc->llx >= 1 && param_loc->lly >= 1 &&
574 param_loc->urx >= 1 && param_loc->ury >= 1, CPL_ERROR_ILLEGAL_INPUT,
575 return CPL_ERROR_ILLEGAL_INPUT,
576 "Coordinates must be strictly positive");
578 cpl_error_ensure(param_loc->urx >= param_loc->llx, CPL_ERROR_ILLEGAL_INPUT,
579 return CPL_ERROR_ILLEGAL_INPUT,
580 "urx (%ld) must be larger equal than llx (%ld)",
581 (
long)param_loc->urx, (
long)param_loc->llx);
582 cpl_error_ensure(param_loc->ury >= param_loc->lly, CPL_ERROR_ILLEGAL_INPUT,
583 return CPL_ERROR_ILLEGAL_INPUT,
584 "ury (%ld) must be larger equal than lly (%ld)",
585 (
long)param_loc->ury, (
long)param_loc->lly);
587 cpl_error_ensure(param_loc->urx <= max_x, CPL_ERROR_ILLEGAL_INPUT,
588 return CPL_ERROR_ILLEGAL_INPUT,
589 "urx %zu larger than maximum %zu",
590 (
size_t)param_loc->urx, (
size_t)max_x);
592 cpl_error_ensure(param_loc->ury <= max_y, CPL_ERROR_ILLEGAL_INPUT,
593 return CPL_ERROR_ILLEGAL_INPUT,
594 "ury %zu larger than maximum %zu",
595 (
size_t)param_loc->ury, (
size_t)max_y);
596 return CPL_ERROR_NONE;
616 const char *base_context,
618 const char *name_prefix,
619 const hdrl_parameter *defaults)
621 cpl_ensure(base_context && prefix && name_prefix && defaults,
622 CPL_ERROR_NULL_INPUT, NULL);
625 CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
627 cpl_parameterlist *parlist = cpl_parameterlist_new();
630 hdrl_setup_vparameter(parlist, prefix,
".", name_prefix,
632 "Lower left x pos. (FITS) defining the region",
635 hdrl_setup_vparameter(parlist, prefix,
".", name_prefix,
637 "Lower left y pos. (FITS) defining the region",
641 hdrl_setup_vparameter(parlist, prefix,
".", name_prefix,
643 "Upper right x pos. (FITS) defining the region",
647 hdrl_setup_vparameter(parlist, prefix,
".", name_prefix,
649 "Upper right y pos. (FITS) defining the region",
652 if (cpl_error_get_code()) {
653 cpl_parameterlist_delete(parlist);
675 const cpl_parameterlist * parlist,
676 const char * base_context,
677 const char * name_prefix)
679 cpl_size llx, lly, urx, ury ;
680 cpl_error_ensure(base_context && parlist, CPL_ERROR_NULL_INPUT,
681 return NULL,
"NULL Input Parameters");
682 const char * sep = strlen(base_context) > 0 ?
"." :
"";
683 const char * points[] = {
"llx",
"lly",
"urx",
"ury"};
684 cpl_size * dest[] = {&llx, &lly, &urx, &ury};
685 for (
size_t i = 0; i < 4; i++) {
686 char * name = cpl_sprintf(
"%s%s%s%s", base_context, sep,
687 name_prefix, points[i]);
688 const cpl_parameter * par = cpl_parameterlist_find_const(parlist, name);
689 *(dest[i]) = cpl_parameter_get_int(par);
693 if (cpl_error_get_code()) {
694 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
695 "Error while parsing parameterlist "
696 "with base_context %s", base_context);
716 hdrl_parameter * rect_region,
720 hdrl_rect_region_parameter * rr_loc =
721 (hdrl_rect_region_parameter *)rect_region ;
723 cpl_error_ensure(rect_region != 0, CPL_ERROR_NULL_INPUT,
724 return CPL_ERROR_NULL_INPUT,
"region input must not be NULL");
726 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
727 "Expected Rect Region parameter") ;
729 if (nx > 0 && rr_loc->llx < 1) rr_loc->llx += nx;
730 if (ny > 0 && rr_loc->lly < 1) rr_loc->lly += ny;
731 if (nx > 0 && rr_loc->urx < 1) rr_loc->urx += nx;
732 if (ny > 0 && rr_loc->ury < 1) rr_loc->ury += ny;
747 if(x == NULL || l <= 1)
return CPL_TRUE;
749 for(cpl_size i = 0; i < l - 1; i++){
750 if(x[i] >= x[i + 1])
return CPL_FALSE;
768 const cpl_size sample_len,
const cpl_boolean sort_decreasing){
770 cpl_propertylist * col_list = cpl_propertylist_new();
772 cpl_propertylist_append_bool(col_list,
"x", sort_decreasing);
774 cpl_table * tb = cpl_table_new(sample_len);
776 cpl_table_wrap_double(tb, x,
"x");
779 cpl_table_wrap_double(tb, y1,
"y1");
782 cpl_table_wrap_double(tb, y2,
"y2");
784 cpl_table_sort(tb, col_list);
786 cpl_table_unwrap(tb,
"x");
789 cpl_table_unwrap(tb,
"y1");
792 cpl_table_unwrap(tb,
"y2");
794 cpl_table_delete(tb);
795 cpl_propertylist_delete(col_list);
816 const char * sep = sep_ ? sep_ :
"";
817 cpl_ensure(n > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
820 for (
int i = 0; i < n; i++) {
822 char * val = va_arg(vl,
char *);
823 if (val == NULL || strlen(val) == 0) {
827 res = cpl_strdup(val);
830 res = cpl_sprintf(
"%s%s%s", res, sep, val);
853int hdrl_get_tempfile(
const char * dir, cpl_boolean unlink)
856 const char * tmpdirs[] = {
864 const size_t ndirs =
sizeof(tmpdirs) /
sizeof(tmpdirs[0]);
865 const char * tmpdir = NULL;
867 if (dir && access(dir, W_OK) == 0) {
871 for (
size_t i = 0; i < ndirs; i++) {
872 if (tmpdirs[i] && access(tmpdirs[i], W_OK) == 0) {
882 int fd = mkstemp(
template);
884 cpl_error_set_message(cpl_func, CPL_ERROR_FILE_IO,
885 "Temporary file creation failed: %s",
891 cpl_msg_debug(cpl_func,
"Created tempfile %s",
template);
909char * hdrl_get_cwd(
void)
917 if (getcwd(buf, n) != 0) {
920 else if (errno == ERANGE) {
928 cpl_error_set_message(cpl_func, CPL_ERROR_FILE_IO,
929 "Could not determine current working "
930 "directory: %s", strerror(errno));
970hdrl_normalize_imagelist_by_vector(
const cpl_vector * scale,
971 const cpl_vector * scale_e,
972 const hdrl_scale_type scale_type,
973 cpl_imagelist * data,
974 cpl_imagelist * errors)
976 cpl_ensure_code(scale, CPL_ERROR_NULL_INPUT);
977 cpl_ensure_code(scale_e, CPL_ERROR_NULL_INPUT);
978 cpl_ensure_code(data, CPL_ERROR_NULL_INPUT);
979 cpl_ensure_code(errors, CPL_ERROR_NULL_INPUT);
980 cpl_ensure_code(cpl_vector_get_size(scale) == cpl_imagelist_get_size(data),
981 CPL_ERROR_ILLEGAL_INPUT);
982 cpl_ensure_code(cpl_vector_get_size(scale_e) == cpl_vector_get_size(scale),
983 CPL_ERROR_ILLEGAL_INPUT);
984 cpl_ensure_code(cpl_imagelist_get_size(errors) ==
985 cpl_imagelist_get_size(data), CPL_ERROR_ILLEGAL_INPUT);
987 for (
size_t i = 1; i < (size_t)cpl_imagelist_get_size(data); i++) {
988 const hdrl_data_t dfirst = cpl_vector_get(scale, 0);
989 const hdrl_error_t efirst = cpl_vector_get(scale_e, 0);
990 cpl_image * dimg = cpl_imagelist_get(data, i);
991 cpl_image * eimg = cpl_imagelist_get(errors, i);
992 if (scale_type == HDRL_SCALE_ADDITIVE) {
993 hdrl_data_t dscale_o = cpl_vector_get(scale, i);
994 hdrl_error_t escale_o = cpl_vector_get(scale_e, i);
997 hdrl_data_t dscale = dfirst;
998 hdrl_error_t escale = efirst;
999 hdrl_elemop_sub(&dscale, &escale, 1, &dscale_o, &escale_o, 1, NULL);
1002 hdrl_elemop_image_add_scalar(dimg, eimg, dscale, escale);
1004 else if (scale_type == HDRL_SCALE_MULTIPLICATIVE) {
1005 hdrl_data_t dscale_o = cpl_vector_get(scale, i);
1006 hdrl_error_t escale_o = cpl_vector_get(scale_e, i);
1008 if (dscale_o == 0.) {
1009 cpl_msg_warning(cpl_func,
"scale factor of image %zu is "
1011 cpl_image_add_scalar(dimg, NAN);
1012 cpl_image_add_scalar(eimg, NAN);
1013 cpl_image_reject_value(dimg, CPL_VALUE_NAN);
1014 cpl_image_reject_value(eimg, CPL_VALUE_NAN);
1019 hdrl_data_t dscale = dfirst;
1020 hdrl_error_t escale = efirst;
1021 hdrl_elemop_div(&dscale, &escale, 1, &dscale_o, &escale_o, 1, NULL);
1024 hdrl_elemop_image_mul_scalar(dimg, eimg, dscale, escale);
1027 return cpl_error_set_message(cpl_func, CPL_ERROR_UNSUPPORTED_MODE,
1028 "Unsupported scale type");
1031 if (cpl_error_get_code())
1035 return cpl_error_get_code();
1071hdrl_normalize_imagelist_by_imagelist(
const cpl_imagelist * scale,
1072 const cpl_imagelist * scale_e,
1073 const hdrl_scale_type scale_type,
1074 cpl_imagelist * data,
1075 cpl_imagelist * errors)
1077 cpl_ensure_code(scale, CPL_ERROR_NULL_INPUT);
1078 cpl_ensure_code(scale_e, CPL_ERROR_NULL_INPUT);
1079 cpl_ensure_code(data, CPL_ERROR_NULL_INPUT);
1080 cpl_ensure_code(errors, CPL_ERROR_NULL_INPUT);
1081 cpl_ensure_code(cpl_imagelist_get_size(scale) ==
1082 cpl_imagelist_get_size(data), CPL_ERROR_ILLEGAL_INPUT);
1083 cpl_ensure_code(cpl_imagelist_get_size(scale_e) ==
1084 cpl_imagelist_get_size(scale), CPL_ERROR_ILLEGAL_INPUT);
1085 cpl_ensure_code(cpl_imagelist_get_size(errors) ==
1086 cpl_imagelist_get_size(data), CPL_ERROR_ILLEGAL_INPUT);
1088 for (
size_t i = 1; i < (size_t)cpl_imagelist_get_size(data); i++) {
1089 cpl_image * dscale =
1090 cpl_image_duplicate(cpl_imagelist_get_const(scale, 0));
1091 cpl_image * escale =
1092 cpl_image_duplicate(cpl_imagelist_get_const(scale_e, 0));
1093 cpl_image * dimg = cpl_imagelist_get(data, i);
1094 cpl_image * eimg = cpl_imagelist_get(errors, i);
1095 const cpl_image * dscale_o = cpl_imagelist_get_const(scale, i);
1096 const cpl_image * escale_o = cpl_imagelist_get_const(scale_e, i);
1098 if (scale_type == HDRL_SCALE_ADDITIVE) {
1100 hdrl_elemop_image_sub_image(dscale, escale, dscale_o, escale_o);
1103 hdrl_elemop_image_add_image(dimg, eimg, dscale, escale);
1105 else if (scale_type == HDRL_SCALE_MULTIPLICATIVE) {
1107 hdrl_elemop_image_div_image(dscale, escale, dscale_o, escale_o);
1110 hdrl_elemop_image_mul_image(dimg, eimg, dscale, escale);
1113 cpl_image_delete(dscale);
1114 cpl_image_delete(escale);
1115 return cpl_error_set_message(cpl_func, CPL_ERROR_UNSUPPORTED_MODE,
1116 "Unsupported scale type");
1119 cpl_image_delete(dscale);
1120 cpl_image_delete(escale);
1122 if (cpl_error_get_code())
1126 return cpl_error_get_code();
1144cpl_vector * hdrl_image_to_vector(
1145 const cpl_image * source,
1146 const cpl_mask * bpm)
1148 cpl_ensure(source != NULL, CPL_ERROR_NULL_INPUT, NULL);
1149 cpl_vector * vec_source = NULL;
1153 CPL_DIAG_PRAGMA_PUSH_IGN(-Wcast-qual);
1154 cpl_image * d_img = cpl_image_get_type(source) == CPL_TYPE_DOUBLE ?
1155 (cpl_image *)source : cpl_image_cast(source, CPL_TYPE_DOUBLE);
1156 CPL_DIAG_PRAGMA_POP;
1158 const cpl_size naxis1 = cpl_image_get_size_x(source);
1159 const cpl_size naxis2 = cpl_image_get_size_y(source);
1162 const double * restrict sdata = cpl_image_get_data_double_const(d_img);
1163 const cpl_binary * restrict bpmd = NULL;
1166 double * restrict ddata = cpl_malloc(naxis1 * naxis2 *
1171 bpmd = cpl_mask_get_data_const(bpm);
1172 else if (cpl_image_get_bpm_const(source) != NULL)
1173 bpmd = cpl_mask_get_data_const(cpl_image_get_bpm_const(source));
1176 memcpy(ddata, sdata, naxis1 * naxis2 *
sizeof(ddata[0]));
1177 j = naxis1 * naxis2;
1181 for (
long i = 0; i < naxis1 * naxis2; i++) {
1182 if (bpmd[i] == CPL_BINARY_0) {
1183 ddata[j] = sdata[i];
1189 assert(j == naxis1 * naxis2 -
1190 (bpm ? cpl_mask_count(bpm) : cpl_image_count_rejected(source)));
1193 vec_source = cpl_vector_wrap(j, ddata);
1199 if (d_img != source) cpl_image_delete(d_img);
1210struct hdrl_vector_cache_ {
1212 cache_bucket buckets[];
1240hdrl_vector_cache * hdrl_vector_cache_new(cpl_size max_cached_size,
1241 cpl_size ncached_entries)
1244 if (max_cached_size > 50)
1246 hdrl_vector_cache * c = cpl_malloc(
sizeof(hdrl_vector_cache) +
1247 sizeof(cache_bucket) * (max_cached_size + 1));
1248 c->nbuckets = max_cached_size + 1;
1249 for (cpl_size i = 0; i < c->nbuckets; i++) {
1250 c->buckets[i].available = 0;
1251 c->buckets[i].nspace = ncached_entries;
1252 c->buckets[i].ptrs = cpl_calloc(
sizeof(cpl_vector*), ncached_entries);
1264void hdrl_vector_cache_delete(hdrl_vector_cache * cache)
1268 for (cpl_size i = 0; i < cache->nbuckets; i++) {
1269 for (
size_t j = 0; j < cache->buckets[i].available; j++) {
1270 cpl_vector_delete(cache->buckets[i].ptrs[j]);
1272 cpl_free(cache->buckets[i].ptrs);
1294cpl_vector * hdrl_cplvector_new_from_cache(hdrl_vector_cache * cache, cpl_size sz)
1297 return cpl_vector_new(sz);
1299 if (sz < cache->nbuckets) {
1300 if (cache->buckets[sz].available > 0) {
1301 return cache->buckets[sz].ptrs[--(cache->buckets[sz].available)];
1304 return cpl_vector_new(sz);
1318void hdrl_cplvector_delete_to_cache(hdrl_vector_cache * cache, cpl_vector * v)
1324 const cpl_size sz = cpl_vector_get_size(v);
1325 if (sz < cache->nbuckets) {
1326 if (cache->buckets[sz].available < cache->buckets[sz].nspace) {
1327 cache->buckets[sz].ptrs[(cache->buckets[sz].available)++] = v;
1332 cpl_vector_delete(v);
1339static cpl_vector * imagelist_to_vector(
const cpl_imagelist * list,
1343 const double ** imgdatabuf,
1344 const cpl_binary ** maskbuf, hdrl_vector_cache * cache)
1346 const long nz = list ? cpl_imagelist_get_size(list) : -1;
1347 unsigned long j = 0;
1348 cpl_vector * vec = hdrl_cplvector_new_from_cache(cache, nz);
1349 double * restrict ddata = cpl_vector_get_data(vec);
1351 if (imgdatabuf && maskbuf) {
1352 for (
long k = 0; k < nz; k++) {
1353 double v = imgdatabuf[k][(y - 1) * nx + (x - 1)];
1354 cpl_binary rej = maskbuf[k] ?
1355 maskbuf[k][(y - 1) * nx + (x - 1)] : 0;
1364 for (
long k = 0; k < nz; k++) {
1366 const cpl_image * img = cpl_imagelist_get_const(list, k);
1367 double v = cpl_image_get(img, x, y, &rej);
1377 if ((cpl_size)j != nz) {
1378 cpl_vector_set_size(vec, j);
1382 hdrl_cplvector_delete_to_cache(cache, vec);
1401cpl_vector * hdrl_imagelist_to_vector(
const cpl_imagelist * list,
1406 const cpl_size nz = list ? cpl_imagelist_get_size(list) : -1;
1407 cpl_ensure(list != NULL, CPL_ERROR_NULL_INPUT, NULL);
1408 cpl_ensure(nz > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
1409 cpl_ensure(x > 0, CPL_ERROR_ACCESS_OUT_OF_RANGE, NULL);
1410 cpl_ensure(y > 0, CPL_ERROR_ACCESS_OUT_OF_RANGE, NULL);
1413 const cpl_image * img = cpl_imagelist_get_const(list, 0);
1414 const cpl_size ny = cpl_image_get_size_y(img);
1415 nx = cpl_image_get_size_x(img);
1416 cpl_ensure(x <= nx, CPL_ERROR_ACCESS_OUT_OF_RANGE, NULL);
1417 cpl_ensure(y <= ny, CPL_ERROR_ACCESS_OUT_OF_RANGE, NULL);
1419 return imagelist_to_vector(list, nx, x, y, NULL, NULL, NULL);
1439cpl_error_code hdrl_imagelist_to_vector_row(
const cpl_imagelist * list,
1442 hdrl_vector_cache * cache)
1446 const long nz = list ? (long)cpl_imagelist_get_size(list) : 0;
1447 cpl_ensure_code(list != NULL, CPL_ERROR_NULL_INPUT);
1448 cpl_ensure_code(nz > 0, CPL_ERROR_ILLEGAL_INPUT);
1449 cpl_ensure_code(y > 0, CPL_ERROR_ACCESS_OUT_OF_RANGE);
1452 const cpl_image * img = cpl_imagelist_get_const(list, 0);
1453 const cpl_size ny = cpl_image_get_size_y(img);
1454 cpl_ensure_code(y <= ny, CPL_ERROR_ACCESS_OUT_OF_RANGE);
1455 nx = (long)cpl_image_get_size_x(img);
1456 isdouble = cpl_image_get_type(img) == CPL_TYPE_DOUBLE;
1459 const double * imgdatabuf[nz];
1460 const cpl_binary * maskbuf[nz];
1461 for (
long i = 0; i < nz && isdouble; i++) {
1462 const cpl_image * img = cpl_imagelist_get_const(list, i);
1463 const cpl_mask * bpm = cpl_image_get_bpm_const(img);
1464 imgdatabuf[i] = cpl_image_get_data_double_const(img);
1466 maskbuf[i] = cpl_mask_get_data_const(bpm);
1472 for (
long x = 0; x < nx; x++) {
1474 out[x] = imagelist_to_vector(list, nx, x + 1, y,
1475 imgdatabuf, maskbuf, cache);
1478 out[x] = imagelist_to_vector(list, nx, x + 1, y, NULL, NULL, cache);
1481 return cpl_error_get_code();
1498hdrl_imagelist_to_cplwrap(
const hdrl_imagelist * list,
1499 cpl_imagelist ** data,
1500 cpl_imagelist ** errs)
1502 cpl_ensure_code(list, CPL_ERROR_NULL_INPUT);
1504 *data = cpl_imagelist_new();
1507 *errs = cpl_imagelist_new();
1519 if (cpl_error_get_code()) {
1521 cpl_imagelist_unwrap(*data);
1525 cpl_imagelist_unwrap(*errs);
1529 return cpl_error_get_code();
1546hdrl_medianfilter_image_grid(
const cpl_image * ima, cpl_matrix * x, cpl_matrix * y,
1547 cpl_size filtersize_x, cpl_size filtersize_y)
1549 cpl_error_ensure(ima != NULL, CPL_ERROR_NULL_INPUT,
return NULL,
1550 "NULL input image");
1551 cpl_error_ensure(filtersize_x > 0 && filtersize_y > 0 ,
1552 CPL_ERROR_INCOMPATIBLE_INPUT,
return NULL,
1553 "All function parameters must be greater then Zero");
1555 const cpl_size nx = cpl_image_get_size_x(ima);
1556 const cpl_size ny = cpl_image_get_size_y(ima);
1557 const cpl_size steps_x = cpl_matrix_get_nrow(x);
1558 const cpl_size steps_y = cpl_matrix_get_nrow(y);
1560 cpl_image * ima_local = cpl_image_new(steps_x, steps_y, CPL_TYPE_DOUBLE);
1562 for (cpl_size iy = 0; iy < steps_y; iy++) {
1563 cpl_size middlep_y = cpl_matrix_get(y, iy, 0);
1564 for (cpl_size ix = 0; ix < steps_x; ix++) {
1565 cpl_size middlep_x = cpl_matrix_get(x, ix, 0);
1567 cpl_size lowerlimit_x = CX_MAX(middlep_x - filtersize_x, 1);
1568 cpl_size lowerlimit_y = CX_MAX(middlep_y - filtersize_y, 1);
1569 cpl_size upperlimit_x = CX_MIN(middlep_x + filtersize_x, nx);
1570 cpl_size upperlimit_y = CX_MIN(middlep_y + filtersize_y, ny);
1572 double median = cpl_image_get_median_window(ima, lowerlimit_x,
1573 lowerlimit_y, upperlimit_x, upperlimit_y);
1575 cpl_image_set(ima_local, ix + 1, iy + 1, median);
1577 cpl_msg_debug(cpl_func,
"middlep_x: %lld, middlep_y: %lld, median: "
1578 "%g", middlep_x, middlep_y, median);
1595cpl_matrix * hdrl_matrix_linspace(
1600 cpl_matrix * x = cpl_matrix_new(stop / step, 1);
1601 for (intptr_t i = 0; start + i * step < stop && i < stop / step; i++) {
1602 cpl_matrix_set(x, i, 0, start + i * step);
1624cpl_matrix * hdrl_fit_legendre(
1628 cpl_matrix * grid_x,
1629 cpl_matrix * grid_y,
1633 cpl_size nx2 = cpl_matrix_get_nrow(grid_x);
1634 cpl_size ny2 = cpl_matrix_get_nrow(grid_y);
1635 cpl_matrix * xpolys =
1637 cpl_matrix * ypolys =
1639 cpl_matrix * tensors =
1642 cpl_matrix * mimage = cpl_matrix_wrap(nx2 * ny2, 1, cpl_image_get_data(img));
1643 cpl_matrix * coeffs = cpl_matrix_solve_normal(tensors, mimage);
1644 cpl_matrix_unwrap(mimage);
1645 cpl_matrix_delete(xpolys);
1646 cpl_matrix_delete(ypolys);
1647 cpl_matrix_delete(tensors);
1665cpl_image * hdrl_legendre_to_image(
1666 cpl_matrix * coeffs,
1674 cpl_matrix * x = hdrl_matrix_linspace(0, nx, 1);
1675 cpl_matrix * y = hdrl_matrix_linspace(0, ny, 1);
1676 cpl_matrix * xpolys =
1678 cpl_matrix * ypolys =
1680 cpl_matrix * tensors =
1683 cpl_matrix * result = cpl_matrix_product_create(tensors, coeffs);
1684 cpl_image * iresult = cpl_image_wrap(nx, ny, CPL_TYPE_DOUBLE,
1685 cpl_matrix_get_data(result));
1686 cpl_matrix_delete(x);
1687 cpl_matrix_delete(y);
1688 cpl_matrix_delete(xpolys);
1689 cpl_matrix_delete(ypolys);
1690 cpl_matrix_delete(tensors);
1691 cpl_matrix_unwrap(result);
1704int hdrl_check_maskequality(
const cpl_mask * mask1,
const cpl_mask * mask2)
1706 cpl_ensure(mask1, CPL_ERROR_NULL_INPUT, 1);
1707 cpl_ensure(mask2, CPL_ERROR_NULL_INPUT, 1);
1709 cpl_size m1nx = cpl_mask_get_size_x(mask1);
1710 cpl_size m1ny = cpl_mask_get_size_y(mask1);
1711 cpl_size m2nx = cpl_mask_get_size_x(mask2);
1712 cpl_size m2ny = cpl_mask_get_size_y(mask2);
1713 cpl_ensure(m1nx == m2nx, CPL_ERROR_NONE, 1);
1714 cpl_ensure(m1ny == m2ny, CPL_ERROR_NONE, 1);
1716 const cpl_binary * dm1bpm = cpl_mask_get_data_const(mask1);
1717 const cpl_binary * dm2bpm = cpl_mask_get_data_const(mask2);
1718 if (memcmp(dm1bpm, dm2bpm, m1nx * m1ny) != 0) {
1726static const cpl_image *
1727image_const_row_view_create(
const cpl_image * img,
1731 const size_t dsz = cpl_type_get_sizeof(cpl_image_get_type(img));
1732 const cpl_size nx = cpl_image_get_size_x(img);
1733 const char * d = cpl_image_get_data_const(img);
1734 size_t offset = (ly - 1) * nx;
1735 cpl_size nny = uy - ly + 1;
1737 CPL_DIAG_PRAGMA_PUSH_IGN(-Wcast-qual);
1739 cpl_image * wimg = cpl_image_wrap(nx, nny, cpl_image_get_type(img),
1740 (
char*)d + offset * dsz);
1742 const cpl_mask * omask = cpl_image_get_bpm_const(img);
1744 cpl_mask * mask = cpl_mask_wrap(nx, nny,
1745 (cpl_binary*)cpl_mask_get_data_const(omask) + offset);
1746 cpl_mask_delete(hcpl_image_set_bpm(wimg, mask));
1749 CPL_DIAG_PRAGMA_POP;
1755image_const_row_view_delete(
const cpl_image * img)
1757 CPL_DIAG_PRAGMA_PUSH_IGN(-Wcast-qual);
1759 cpl_mask_unwrap(cpl_image_unset_bpm((cpl_image*)(img)));
1760 cpl_image_unwrap((cpl_image *)img);
1762 CPL_DIAG_PRAGMA_POP;
1783hdrl_parallel_filter_image(
const cpl_image * img,
1784 const cpl_matrix * kernel,
1785 const cpl_mask * mask,
1786 const cpl_filter_mode mode)
1788 cpl_ensure(img, CPL_ERROR_NULL_INPUT, NULL);
1789 intptr_t nx = cpl_image_get_size_x(img);
1790 intptr_t ny = cpl_image_get_size_y(img);
1793 const cpl_border_mode border = CPL_BORDER_FILTER;
1794 cpl_ensure((kernel && !mask) || (!kernel && mask),
1795 CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
1797 ky = cpl_matrix_get_nrow(kernel);
1798 kx = cpl_matrix_get_ncol(kernel);
1801 ky = cpl_mask_get_size_y(mask);
1802 kx = cpl_mask_get_size_x(mask);
1804 cpl_ensure(ky % 2 == 1, CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
1805 cpl_ensure(ny >= ky, CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
1806 cpl_ensure(nx >= kx, CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
1808 intptr_t hk = (ky / 2);
1809 cpl_image * res = cpl_image_new(nx, ny, cpl_image_get_type(img));
1811 cpl_image_get_bpm(res);
1815 const cpl_image * slice = image_const_row_view_create(img, 1, ky);
1816 cpl_image * slres = cpl_image_duplicate(slice);
1818 cpl_image_filter(slres, slice, kernel, mode, border);
1821 cpl_image_filter_mask(slres, slice, mask, mode, border);
1824 const cpl_image * slice2 = image_const_row_view_create(slres, 1, hk);
1825 cpl_image_copy(res, slice2, 1, 1);
1826 image_const_row_view_delete(slice2);
1827 image_const_row_view_delete(slice);
1828 cpl_image_delete(slres);
1831 const intptr_t s = 200;
1833HDRL_OMP(omp parallel
for lastprivate(y)
if (ny > s + ky))
1834 for (y = hk; y < ny - ky - (ny - ky) % s; y+=s) {
1835 intptr_t l = (y + 1) - hk;
1836 intptr_t u = (y + 1 + s) + hk - 1;
1837 const cpl_image * slice = image_const_row_view_create(img, l, u);
1838 cpl_image * slres = cpl_image_new(nx, u - l + 1,
1839 cpl_image_get_type(slice));
1841 cpl_image_filter(slres, slice, kernel, mode, border);
1844 cpl_image_filter_mask(slres, slice, mask, mode, border);
1846 const cpl_image * slice2 =
1847 image_const_row_view_create(slres, hk + 1, hk + s);
1848 cpl_image_copy(res, slice2, 1, y + 1);
1849 image_const_row_view_delete(slice);
1850 image_const_row_view_delete(slice2);
1851 cpl_image_delete(slres);
1854 if (y + 1 - hk < ny) {
1855 const cpl_image * slice = image_const_row_view_create(img, y + 1 - hk, ny);
1856 cpl_image * slres = cpl_image_duplicate(slice);
1858 cpl_image_filter(slres, slice, kernel, mode, border);
1861 cpl_image_filter_mask(slres, slice, mask, mode, border);
1863 const cpl_image * slice2 =
1864 image_const_row_view_create(slres, hk + 1, cpl_image_get_size_y(slice));
1865 cpl_image_copy(res, slice2, 1, y + 1);
1866 image_const_row_view_delete(slice);
1867 image_const_row_view_delete(slice2);
1868 cpl_image_delete(slres);
1893hdrl_wcs_convert(
const cpl_wcs *wcs,
const cpl_matrix *from,
1894 cpl_matrix **to, cpl_array **status,
1895 cpl_wcs_trans_mode transform)
1897 size_t nr = cpl_matrix_get_nrow(from);
1898 size_t nc = cpl_matrix_get_ncol(from);
1899 const size_t s = 4000;
1902 cpl_ensure_code(to, CPL_ERROR_NULL_INPUT);
1903 cpl_ensure_code(status, CPL_ERROR_NULL_INPUT);
1904 cpl_ensure_code(wcs, CPL_ERROR_NULL_INPUT);
1905 cpl_ensure_code(from, CPL_ERROR_NULL_INPUT);
1907 *status = cpl_array_new(nr, CPL_TYPE_INT);
1908 cpl_ensure_code(*status, CPL_ERROR_NULL_INPUT);
1909 dstatus = cpl_array_get_data_int(*status);
1910 *to = cpl_matrix_new(nr, nc);
1911 cpl_error_code err = CPL_ERROR_NONE;
1913HDRL_OMP(omp parallel
for if (nr > s))
1914 for (i = 0; i < nr; i+=s) {
1915 cpl_matrix * lfrom = cpl_matrix_extract(from, i, 0, 1, 1, CPL_MIN(s, nr - i), nc);
1916 cpl_matrix * lto = NULL;
1917 cpl_array * lstatus = NULL;
1918 cpl_error_code lerr = cpl_wcs_convert(wcs, lfrom, <o, &lstatus, transform);
1920 cpl_matrix_copy(*to, lto, i, 0);
1922 if (lstatus != NULL) {
1923 memcpy(&dstatus[i], cpl_array_get_data_int(lstatus),
1924 cpl_array_get_size(lstatus) *
sizeof(*dstatus));
1926 cpl_array_delete(lstatus);
1927 cpl_matrix_delete(lfrom);
1928 cpl_matrix_delete(lto);
1929 if (lerr != CPL_ERROR_NONE) {
1930HDRL_OMP(omp critical(hdrl_hdrlwcserror))
1935 if (err == CPL_ERROR_UNSUPPORTED_MODE) {
1936 cpl_matrix_delete(*to);
1938 cpl_array_delete(*status);
1942 return cpl_error_set(cpl_func, err);
1957cpl_mask * hcpl_image_set_bpm(cpl_image * self, cpl_mask * bpm)
1959#if CPL_VERSION_CODE < CPL_VERSION(6, 4, 0)
1962 cpl_ensure(self != NULL, CPL_ERROR_NULL_INPUT, NULL);
1969 return cpl_image_set_bpm(self, bpm);
1973double hcpl_vector_get_mad_window(cpl_vector * vec,
1979 struct _cpl_image_ img;
1980 img.pixels = cpl_vector_get_data(vec);
1981 img.nx = cpl_vector_get_size(vec);
1984 img.type = CPL_TYPE_DOUBLE;
1985 return cpl_image_get_mad_window(&img, llx, 1, urx, 1, sigma);
1988double hcpl_gaussian_eval_2d(
const cpl_array * self,
double x,
double y)
1990#if CPL_VERSION_CODE < CPL_VERSION(6, 4, 0)
1991 cpl_errorstate prestate = cpl_errorstate_get();
1992 const double B = cpl_array_get_double(self, 0, NULL);
1993 const double A = cpl_array_get_double(self, 1, NULL);
1994 const double R = cpl_array_get_double(self, 2, NULL);
1995 const double M_x = cpl_array_get_double(self, 3, NULL);
1996 const double M_y = cpl_array_get_double(self, 4, NULL);
1997 const double S_x = cpl_array_get_double(self, 5, NULL);
1998 const double S_y = cpl_array_get_double(self, 6, NULL);
2002 if (!cpl_errorstate_is_equal(prestate)) {
2003 (void)cpl_error_set_where(cpl_func);
2004 }
else if (cpl_array_get_size(self) != 7) {
2005 (void)cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_INPUT);
2006 }
else if (fabs(R) < 1.0 && S_x != 0.0 && S_y != 0.0) {
2007 const double x_n = (x - M_x) / S_x;
2008 const double y_n = (y - M_y) / S_y;
2010 value = B + A / (CPL_MATH_2PI * S_x * S_y * sqrt(1 - R * R)) *
2011 exp(-0.5 / (1 - R * R) * ( x_n * x_n + y_n * y_n
2012 - 2.0 * R * x_n * y_n));
2013 }
else if (fabs(R) > 1.0) {
2014 (void)cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_OUTPUT,
2015 "fabs(R=%g) > 1", R);
2017 (void)cpl_error_set_message(cpl_func, CPL_ERROR_DIVISION_BY_ZERO,
2018 "R=%g. Sigma=(%g, %g)", R, S_x, S_y);
2023 return cpl_gaussian_eval_2d(self, x, y);
2070 hdrl_value aRA, hdrl_value aDEC, hdrl_value aLST,
2071 hdrl_value aExptime, hdrl_value aLatitude,
2072 hdrl_airmass_approx type)
2075 hdrl_value retErr = {-1.,0.};
2076 cpl_ensure( aRA.data >= 0. && aRA.data < 360. && aRA.error >= 0.
2077 && aDEC.data >= -90. && aDEC.data <= 90. && aDEC.error >= 0.
2078 && aLST.data >= 0. && aLST.data < 86400. && aLST.error >= 0.
2079 && aExptime.data >= 0. && aExptime.error >= 0.
2080 && aLatitude.data >= -90. && aLatitude.data <= 90. && aLatitude.error >= 0.
2081 && (type == HDRL_AIRMASS_APPROX_HARDIE || type == HDRL_AIRMASS_APPROX_YOUNG_IRVINE || type == HDRL_AIRMASS_APPROX_YOUNG),
2082 CPL_ERROR_ILLEGAL_INPUT,
2086 hdrl_value HA = {aLST.data * 15. / 3600. - aRA.data,
2087 aLST.error * fabs(15. / 3600.) + aRA.error * fabs(-1.)};
2090 if (HA.data < -180.) HA.data += 360.;
2091 if (HA.data > 180.) HA.data -= 360.;
2095 hdrl_value delta = {aDEC.data * CPL_MATH_RAD_DEG,
2096 aDEC.error * fabs(CPL_MATH_RAD_DEG)};
2098 hdrl_value latitude = {aLatitude.data * CPL_MATH_RAD_DEG,
2099 aLatitude.error * fabs(CPL_MATH_RAD_DEG)};
2101 hdrl_value hourangle = {HA.data * CPL_MATH_RAD_DEG,
2102 HA.error * fabs(CPL_MATH_RAD_DEG)};
2109 hdrl_value cosz = hdrl_get_zenith_distance(hourangle, delta, latitude);
2110 hdrl_value z = {0., cosz.error};
2111 double zlimit = 80.;
2113 if (type == HDRL_AIRMASS_APPROX_HARDIE) {
2115 z.data = acos(cosz.data) * CPL_MATH_DEG_RAD;
2116 z.error = cosz.error * fabs(-CPL_MATH_DEG_RAD / sqrt(1. + pow(cosz.data, 2)));
2118 cpl_error_ensure(z.data <= zlimit, CPL_ERROR_ILLEGAL_OUTPUT,
return retErr,
2119 "Zenith angle %f+-[%f] > %f!", z.data, z.error, zlimit);
2123 cosz.data != 0. && fabs(1./cosz.data) >= FLT_EPSILON && acos(cosz.data) <= CPL_MATH_PI_2,
2124 CPL_ERROR_ILLEGAL_OUTPUT,
2126 "Airmass computation unsuccessful. Object is below the horizon at start (z = %f). Return the error",
2127 cosz.error * fabs(-CPL_MATH_DEG_RAD / sqrt(1. + pow(cosz.data, 2))) );
2129 hdrl_value airmass = {0.,0.};
2130 hdrl_value secansZdist = {1. / cosz.data, cosz.error * fabs(-1. / pow(cosz.data, 2))};
2133 case HDRL_AIRMASS_APPROX_HARDIE:
2134 airmass = hdrl_get_airmass_hardie(secansZdist);
2136 case HDRL_AIRMASS_APPROX_YOUNG_IRVINE:
2137 airmass = hdrl_get_airmass_youngirvine(secansZdist);
2139 case HDRL_AIRMASS_APPROX_YOUNG:
2140 airmass = hdrl_get_airmass_young(cosz);
2147 if (aExptime.data > 0.) {
2149 const double weights[] = {1./6., 2./3., 1./6.};
2150 const cpl_size nweights =
sizeof(weights) /
sizeof(
double);
2152 hdrl_value timeStep = {aExptime.data / (nweights - 1.) * 15. / 3600. * CPL_MATH_RAD_DEG, 0.};
2153 timeStep.error = aExptime.error * fabs( 1. / (nweights - 1.) * 15. / 3600. * CPL_MATH_RAD_DEG);
2155 airmass.data *= weights[0];
2156 airmass.error *= fabs(weights[0]);
2158 for (cpl_size i = 1; i < nweights; i++) {
2160 hdrl_value aux_hourangle = { hourangle.data + i * timeStep.data,
2161 hourangle.error + i * timeStep.error };
2163 cosz = hdrl_get_zenith_distance(aux_hourangle, delta, latitude);
2165 if (type == HDRL_AIRMASS_APPROX_HARDIE) {
2167 z = (hdrl_value){acos(cosz.data) * CPL_MATH_DEG_RAD, 0.};
2168 z.error = cosz.error * fabs(-CPL_MATH_DEG_RAD / sqrt(1. + pow(cosz.data, 2)));
2170 cpl_error_ensure(z.data <= zlimit, CPL_ERROR_ILLEGAL_OUTPUT,
return retErr,
2171 "Zenith angle %f+-[%f] > %f!", z.data, z.error, zlimit);
2175 cosz.data != 0. && fabs(1./cosz.data) >= FLT_EPSILON && acos(cosz.data) <= CPL_MATH_PI_2,
2176 CPL_ERROR_ILLEGAL_OUTPUT,
2178 "timeStep. Object is below the horizon at %s exposure (z=%f).",
2179 i==1 ?
"mid. Return the error" :
"end. Return the error",
2180 cosz.error * fabs(-CPL_MATH_DEG_RAD / sqrt(1. + pow(cosz.data, 2))) );
2182 hdrl_value weight = {0., 0.};
2183 secansZdist = (hdrl_value){1. / cosz.data, cosz.error * fabs(-1. / pow(cosz.data, 2))};
2186 case HDRL_AIRMASS_APPROX_HARDIE:
2187 weight = hdrl_get_airmass_hardie(secansZdist);
2189 case HDRL_AIRMASS_APPROX_YOUNG_IRVINE:
2190 weight = hdrl_get_airmass_youngirvine(secansZdist);
2192 case HDRL_AIRMASS_APPROX_YOUNG:
2193 weight = hdrl_get_airmass_young(cosz);
2196 airmass.data += weights[i] *weight.data;
2197 airmass.error += weights[i] *weight.error;
2201 if (type == HDRL_AIRMASS_APPROX_YOUNG_IRVINE) {
2204 const double airmasslimit = 4.;
2206 cpl_error_ensure(airmass.data <= airmasslimit, CPL_ERROR_ILLEGAL_OUTPUT,
2207 return retErr,
"Airmass larger than %f", airmasslimit);
2230hdrl_value hdrl_get_zenith_distance(
2231 hdrl_value aHourAngle, hdrl_value aDelta, hdrl_value aLatitude)
2234 hdrl_value p0 = {sin(aLatitude.data) * sin(aDelta.data),
2235 aLatitude.error * fabs( cos(aLatitude.data) * sin(aDelta.data)
2236 + aDelta.error * fabs( sin(aLatitude.data) * cos(aDelta.data)))};
2238 hdrl_value p1 = {cos(aLatitude.data) * cos(aDelta.data),
2239 aLatitude.error * fabs(-sin(aLatitude.data) * cos(aDelta.data)
2240 + aDelta.error * fabs(-cos(aLatitude.data) * sin(aDelta.data)))};
2242 hdrl_value z = {p0.data + cos(aHourAngle.data) * p1.data,
2244 + aHourAngle.error * fabs(-sin(aHourAngle.data) * p1.data)
2245 + p1.error * fabs( cos(aHourAngle.data))};
2247 return fabs(z.data) < FLT_EPSILON ? (hdrl_value){0.,0.} : z;
2264hdrl_value hdrl_get_airmass_hardie(hdrl_value hvaSecZ)
2266 double aSecZ = hvaSecZ.data,
2267 aSecZErr = hvaSecZ.error;
2269 hdrl_value secm1 = {aSecZ - 1, aSecZErr};
2271 hdrl_value airmass = {aSecZ - 0.0018167 * secm1.data
2272 - 0.002875 * secm1.data * secm1.data
2273 - 0.0008083 * secm1.data * secm1.data * secm1.data,
2274 aSecZErr +secm1.error * fabs( - 0.0018167
2275 - 2. * 0.002875 * secm1.data
2276 - 3. * 0.0008083 * secm1.data * secm1.data)};
2294hdrl_value hdrl_get_airmass_youngirvine(hdrl_value hvaSecZ)
2296 double aSecZ = hvaSecZ.data,
2297 aSecZErr = hvaSecZ.error;
2299 hdrl_value airmass = {aSecZ * (1. - 0.0012 * (pow(aSecZ, 2) - 1.)),
2300 aSecZErr * fabs( (1. - 0.0012 * (pow(aSecZ, 2) - 1.)) - 2. * 0.0012 * pow(aSecZ, 2))};
2320hdrl_value hdrl_get_airmass_young(hdrl_value hvaCosZt)
2322 double aCosZt = hvaCosZt.data,
2323 aCosZtErr = hvaCosZt.error;
2325 hdrl_value airmass = {(1.002432 * aCosZt * aCosZt + 0.148386 * aCosZt + 0.0096467)
2326 / (aCosZt * aCosZt * aCosZt + 0.149864 * aCosZt * aCosZt + 0.0102963 * aCosZt
2328 aCosZtErr * fabs( ( (2. * 1.002432 * aCosZt + 0.148386) * (aCosZt * aCosZt * aCosZt + 0.149864 * aCosZt * aCosZt + 0.0102963 * aCosZt + 0.000303978)
2329 - (3. * aCosZt * aCosZt + 2. * 0.149864 * aCosZt + 0.0102963) * (1.002432 * aCosZt * aCosZt + 0.148386 * aCosZt + 0.0096467)
2330 ) / pow(aCosZt * aCosZt * aCosZt + 0.149864 * aCosZt * aCosZt + 0.0102963 * aCosZt + 0.000303978, 2) )};
cpl_image * hdrl_image_get_error(hdrl_image *himg)
get error as cpl image
cpl_image * hdrl_image_get_image(hdrl_image *himg)
get data as cpl image
cpl_size hdrl_imagelist_get_size(const hdrl_imagelist *himlist)
Get the number of images in the imagelist.
hdrl_image * hdrl_imagelist_get(const hdrl_imagelist *himlist, cpl_size inum)
Get an image from a list of images.
cpl_matrix * hdrl_mime_linalg_pairwise_column_tensor_products_create(const cpl_matrix *mat1, const cpl_matrix *mat2)
Create selected pairwise tensor products of the columns of two matrices.
cpl_matrix * hdrl_mime_legendre_polynomials_create(int npoly, double a, double b, const cpl_matrix *x)
Create the Legendre polynomial basis on the interval (a,b).
const char * hdrl_get_license(void)
Get the pipeline copyright and license.
void hdrl_sort_on_x(double *x, double *y1, double *y2, const cpl_size sample_len, const cpl_boolean sort_decreasing)
sort in increasing or decreasing order x. Keep aligned with y1 and y2.
hdrl_value hdrl_utils_airmass(hdrl_value aRA, hdrl_value aDEC, hdrl_value aLST, hdrl_value aExptime, hdrl_value aLatitude, hdrl_airmass_approx type)
Compute the effective airmass of an observation. Takes in count the error propagation if you enter th...
cpl_boolean hdrl_is_strictly_monotonic_increasing(const double *x, cpl_size l)
returns CPL_TRUE if x is strictly monotonic increasing
cpl_error_code hdrl_rect_region_parameter_verify(const hdrl_parameter *param, const cpl_size max_x, const cpl_size max_y)
Verify basic correctness of the parameters.
cpl_table * hdrl_eop_data_totable(const char *eop_data, cpl_size data_length)
Export a raw string buffer containing EOP data to a CPL table.
cpl_error_code hdrl_rect_region_parameter_update(hdrl_parameter *rect_region, cpl_size llx, cpl_size lly, cpl_size urx, cpl_size ury)
Update Rect Region Parameters object.
cpl_parameterlist * hdrl_rect_region_parameter_create_parlist(const char *base_context, const char *prefix, const char *name_prefix, const hdrl_parameter *defaults)
Create parameter list for hdrl_rect_region.
cpl_error_code hdrl_rect_region_fix_negatives(hdrl_parameter *rect_region, const cpl_size nx, const cpl_size ny)
wrap negative or zero coordinates around full image size
hdrl_parameter * hdrl_rect_region_parameter_create(cpl_size llx, cpl_size lly, cpl_size urx, cpl_size ury)
Creates Rect Region Parameters object.
cpl_size hdrl_rect_region_get_llx(const hdrl_parameter *p)
get lower left x coordinate of rectangual region
cpl_size hdrl_rect_region_get_urx(const hdrl_parameter *p)
get upper right x coordinate of rectangular region
cpl_size hdrl_rect_region_get_lly(const hdrl_parameter *p)
get lower left y coordinate of rectangual region
char * hdrl_join_string(const char *sep_, int n,...)
join strings together
cpl_size hdrl_rect_region_get_ury(const hdrl_parameter *p)
get upper right y coordinate of rectangual region
hdrl_parameter * hdrl_rect_region_parameter_parse_parlist(const cpl_parameterlist *parlist, const char *base_context, const char *name_prefix)
parse parameterlist for rectangle parameters
cpl_boolean hdrl_rect_region_parameter_check(const hdrl_parameter *self)
Check that the parameter is hdrl_rect_region parameter.