00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022 #ifdef HAVE_CONFIG_H
00023 #include <config.h>
00024 #endif
00025
00026
00027
00028
00029 #include <errno.h>
00030 #include <stdlib.h>
00031 #include <string.h>
00032 #include <strings.h>
00033 #include <math.h>
00034 #include <cpl.h>
00035 #include <muse.h>
00036 #include "muse_exp_align_z.h"
00037
00038
00039 typedef struct muse_indexpair muse_indexpair;
00040
00041 struct muse_indexpair {
00042 cpl_size first;
00043 cpl_size second;
00044 };
00045
00046 typedef void (*muse_free_func)(void *);
00047
00048
00049
00050 static const double deg2as = 3600.;
00051
00052
00053
00054
00055
00056
00057 static void
00058 muse_vfree(void **array, cpl_size size, muse_free_func deallocator)
00059 {
00060 if (array) {
00061 cpl_size idx;
00062 for (idx = 0; idx < size; ++idx) {
00063 if (deallocator) {
00064 deallocator(array[idx]);
00065 }
00066 }
00067 cpl_free(array);
00068 }
00069 return;
00070 }
00071
00072 static int
00073 _muse_condition_less_than(double aValue1, double aValue2)
00074 {
00075 return (aValue1 < aValue2) ? TRUE : FALSE;
00076 }
00077
00078
00092
00093 static cpl_error_code
00094 muse_utils_replace_nan(muse_image *aImage, float aValue)
00095 {
00096 cpl_ensure(aImage, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
00097 cpl_ensure(aImage->data, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_ILLEGAL_INPUT);
00098
00099 cpl_size npixel = cpl_image_get_size_x(aImage->data) *
00100 cpl_image_get_size_y(aImage->data);
00101
00102 float *data = cpl_image_get_data_float(aImage->data);
00103 if (cpl_error_get_code() == CPL_ERROR_TYPE_MISMATCH) {
00104 cpl_error_set_message(__func__, CPL_ERROR_TYPE_MISMATCH,
00105 "MUSE image data buffer has invalid type!");
00106 return CPL_ERROR_TYPE_MISMATCH;
00107 }
00108
00109 cpl_size ipixel;
00110 for (ipixel = 0; ipixel < npixel; ++ipixel) {
00111 if (isnan(data[ipixel])) {
00112 data[ipixel] = aValue;
00113 }
00114 }
00115
00116 if (aImage->stat) {
00117 npixel = cpl_image_get_size_x(aImage->stat) *
00118 cpl_image_get_size_y(aImage->stat);
00119
00120 data = cpl_image_get_data_float(aImage->stat);
00121 if (cpl_error_get_code() == CPL_ERROR_TYPE_MISMATCH) {
00122 cpl_error_set_message(__func__, CPL_ERROR_TYPE_MISMATCH,
00123 "MUSE image data buffer has invalid type!");
00124 return CPL_ERROR_TYPE_MISMATCH;
00125 }
00126
00127 for (ipixel = 0; ipixel < npixel; ++ipixel) {
00128 if (isnan(data[ipixel])) {
00129 data[ipixel] = aValue;
00130 }
00131 }
00132 }
00133
00134 return CPL_ERROR_NONE;
00135 }
00136
00137
00138 #if 0
00139
00155
00156 static muse_image *
00157 muse_fov_load(const char *aFilename)
00158 {
00159 muse_image *image = muse_image_new();
00160
00161 image->header = cpl_propertylist_load(aFilename, 0);
00162 if (!image->header) {
00163 cpl_error_set_message(__func__, cpl_error_get_code(), "Loading primary FITS "
00164 "header of \"%s\" did not succeed", aFilename);
00165 muse_image_delete(image);
00166 return NULL;
00167 }
00168
00169
00170 int extension = cpl_fits_find_extension(aFilename, EXTNAME_DATA);
00171 cpl_propertylist *hdata = cpl_propertylist_load(aFilename, extension);
00172 while (muse_pfits_get_naxis(hdata, 0) != 2) {
00173
00174
00175 cpl_msg_debug(__func__, "Skipping extension %d [%s]", extension,
00176 muse_pfits_get_extname(hdata));
00177 cpl_propertylist_delete(hdata);
00178 extension++;
00179 hdata = cpl_propertylist_load(aFilename, extension);
00180 }
00181 cpl_msg_debug(__func__, "Taking extension %d [%s]", extension,
00182 muse_pfits_get_extname(hdata));
00183
00184
00185 char *extname = cpl_strdup(muse_pfits_get_extname(hdata));
00186 image->data = cpl_image_load(aFilename, CPL_TYPE_FLOAT, 0, extension);
00187 if (!image->data) {
00188 cpl_error_set_message(__func__, MUSE_ERROR_READ_DATA, "Could not load "
00189 "extension %s from \"%s\"", extname, aFilename);
00190 muse_image_delete(image);
00191 cpl_free(extname);
00192 return NULL;
00193 }
00194
00195 if (cpl_propertylist_has(hdata, "BUNIT")) {
00196 cpl_propertylist_append_string(image->header, "BUNIT",
00197 muse_pfits_get_bunit(hdata));
00198 cpl_propertylist_set_comment(image->header, "BUNIT",
00199 cpl_propertylist_get_comment(hdata, "BUNIT"));
00200 } else {
00201 cpl_msg_warning(__func__, "No BUNIT given in extension %d [%s] of \"%s\"!",
00202 extension, extname, aFilename);
00203 }
00204
00205 if (!cpl_propertylist_has(hdata, "CUNIT1") ||
00206 !cpl_propertylist_has(hdata, "CUNIT2")) {
00207 cpl_msg_warning(__func__, "No WCS found in extension %d [%s] of \"%s\"!",
00208 extension, extname, aFilename);
00209 }
00210
00211
00212
00213 cpl_propertylist_erase_regexp(hdata, "(^ESO |" MUSE_WCS_KEYS ")", 1);
00214 cpl_propertylist_append(image->header, hdata);
00215 cpl_propertylist_delete(hdata);
00216
00217
00218
00219 if (extname && !strncmp(extname, EXTNAME_DATA, strlen(EXTNAME_DATA)+1)) {
00220
00221 extension = cpl_fits_find_extension(aFilename, EXTNAME_STAT);
00222 } else {
00223
00224 char *statname = cpl_sprintf("%s_STAT", extname);
00225 extension = cpl_fits_find_extension(aFilename, statname);
00226 cpl_free(statname);
00227 }
00228 if (extension > 0) {
00229 cpl_errorstate status = cpl_errorstate_get();
00230 image->stat = cpl_image_load(aFilename, CPL_TYPE_INT, 0, extension);
00231 if (!cpl_errorstate_is_equal(status)) {
00232 if (cpl_error_get_code() == CPL_ERROR_DATA_NOT_FOUND) {
00233 cpl_errorstate_set(status);
00234 cpl_msg_debug(__func__, "Ignoring empty extension %s in \"%s\"",
00235 EXTNAME_STAT, aFilename);
00236 } else {
00237 cpl_error_set_message(__func__, MUSE_ERROR_READ_STAT, "Could not load "
00238 "extension %s from \"%s\"", EXTNAME_STAT, aFilename);
00239 muse_image_delete(image);
00240 cpl_free(extname);
00241 return NULL;
00242 }
00243 }
00244 }
00245
00246
00247
00248 if (extname && !strncmp(extname, EXTNAME_DATA, strlen(EXTNAME_DATA)+1)) {
00249
00250 extension = cpl_fits_find_extension(aFilename, EXTNAME_DQ);
00251 } else {
00252
00253 char *dqname = cpl_sprintf("%s_DQ", extname);
00254 extension = cpl_fits_find_extension(aFilename, dqname);
00255 cpl_free(dqname);
00256 }
00257 if (extension > 0) {
00258 cpl_errorstate status = cpl_errorstate_get();
00259 image->dq = cpl_image_load(aFilename, CPL_TYPE_INT, 0, extension);
00260 if (!cpl_errorstate_is_equal(status)) {
00261 cpl_errorstate_set(status);
00262 cpl_error_set_message(__func__, MUSE_ERROR_READ_DQ, "Could not load "
00263 "extension %s from \"%s\"", EXTNAME_DQ, aFilename);
00264 muse_image_delete(image);
00265 cpl_free(extname);
00266 return NULL;
00267 }
00268 muse_image_dq_to_nan(image);
00269 }
00270
00271 cpl_free(extname);
00272 return image;
00273 }
00274 #endif
00275
00276
00277
00283
00284 static muse_imagelist *
00285 muse_processing_fov_load_all(muse_processing *aProcessing)
00286 {
00287 cpl_ensure(aProcessing, CPL_ERROR_NULL_INPUT, NULL);
00288 cpl_size nframes = cpl_frameset_get_size(aProcessing->inframes);
00289 cpl_ensure(nframes, CPL_ERROR_DATA_NOT_FOUND, NULL);
00290
00291 muse_imagelist *fovimages = muse_imagelist_new();
00292
00293 cpl_size iexposure = 0;
00294 cpl_size iframe;
00295 for (iframe = 0; iframe < nframes; ++iframe) {
00296 const cpl_frame *frame = cpl_frameset_get_position(aProcessing->inframes,
00297 iframe);
00298 const char *tag = cpl_frame_get_tag(frame);
00299
00300 if (muse_processing_check_intags(aProcessing, tag, strlen(tag))) {
00301 const char *filename = cpl_frame_get_filename(frame);
00302
00303 cpl_msg_debug(__func__, "Loading FOV image '%s' as exposure %"
00304 CPL_SIZE_FORMAT, filename, iexposure + 1);
00305
00306 muse_image *image = muse_fov_load(filename);
00307 if (!image) {
00308 cpl_msg_error(__func__, "Could not load FOV image '%s'", filename);
00309 muse_imagelist_delete(fovimages);
00310 return NULL;
00311 }
00312
00313 muse_imagelist_set(fovimages, image, iexposure);
00314 muse_processing_append_used(aProcessing, (cpl_frame *)frame,
00315 CPL_FRAME_GROUP_RAW, 1);
00316 ++iexposure;
00317 }
00318 }
00319 return fovimages;
00320 }
00321
00322
00329
00330 static cpl_boolean
00331 muse_align_check_detection_params(muse_exp_align_params_t *aParams)
00332 {
00333 if (aParams->threshold <= 0.) {
00334 cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
00335 "Detection threshold must be greater than zero!");
00336 return CPL_FALSE;
00337 }
00338 if (aParams->fwhm < 0.5) {
00339 cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
00340 "FWHM must be greater than 0.5 pixels!");
00341 return CPL_FALSE;
00342 }
00343 if (aParams->srcmax <= aParams->srcmin) {
00344 cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
00345 "Maximum number of sources must be grater than the "
00346 "minimum number of sources!");
00347 return CPL_FALSE;
00348 }
00349 if (aParams->roundmax < aParams->roundmin) {
00350 cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
00351 "Upper limit of the point-source roundness must be "
00352 "greater than the lower limit!");
00353 return CPL_FALSE;
00354 }
00355 if (aParams->sharpmin <= 0.) {
00356 cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
00357 "Low limit of the point-source sharpness must be "
00358 "greater than zero!");
00359 return CPL_FALSE;
00360 }
00361 if (aParams->sharpmax < aParams->sharpmin) {
00362 cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
00363 "Upper limit of the point-source sharpness must be "
00364 "greater than the lower limit!");
00365 return CPL_FALSE;
00366 }
00367 return CPL_TRUE;
00368 }
00369
00370
00378
00379 static cpl_array *
00380 muse_align_parse_valuelist(const char *aValuelist)
00381 {
00382 cpl_ensure(aValuelist, CPL_ERROR_NULL_INPUT, NULL);
00383 if (strlen(aValuelist) == 0) {
00384 cpl_error_set_message(__func__, CPL_ERROR_DATA_NOT_FOUND,
00385 "List of values is empty!");
00386 return NULL;
00387 }
00388
00389 cpl_array *strings = muse_cplarray_new_from_delimited_string(aValuelist, ",");
00390 cpl_size nvalues = cpl_array_get_size(strings);
00391
00392 cpl_array *values = cpl_array_new(nvalues, CPL_TYPE_DOUBLE);
00393 cpl_size ivalue;
00394 for (ivalue = 0; ivalue < nvalues; ++ivalue) {
00395 const char *svalue = cpl_array_get_string(strings, ivalue);
00396 if ((strncasecmp(svalue, "inf", 3) == 0) ||
00397 (strncasecmp(svalue, "nan", 3) == 0)) {
00398 cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
00399 "Illegal value \"%s\" encountered!", svalue);
00400 cpl_array_delete(values);
00401 cpl_array_delete(strings);
00402 return NULL;
00403 }
00404 else {
00405 char *last = NULL;
00406 double value = strtod(svalue, &last);
00407 if( (errno == ERANGE) || ((value == 0.) && (last == svalue))) {
00408 cpl_error_set_message(__func__, CPL_ERROR_TYPE_MISMATCH,
00409 "Could not convert \"%s\" to a decimal number!", svalue);
00410 cpl_array_delete(values);
00411 cpl_array_delete(strings);
00412 return NULL;
00413 }
00414 cpl_array_set_double(values, ivalue, value);
00415 }
00416 }
00417 cpl_array_delete(strings);
00418 return values;
00419 }
00420
00421
00430
00431 static int
00432 muse_align_celestial_from_pixel(cpl_table *aTable,
00433 cpl_propertylist *aWCSHeader)
00434 {
00435 muse_wcs *wcs = muse_wcs_new(aWCSHeader);
00436
00437 const char *unit1 = muse_pfits_get_cunit(aWCSHeader, 1);
00438 const char *unit2 = muse_pfits_get_cunit(aWCSHeader, 2);
00439
00440 if (unit1 && unit2) {
00441 if ((strncmp(unit1, unit2, 4) == 0) && (strncmp(unit1, "deg", 3) == 0)) {
00442 wcs->iscelsph = CPL_TRUE;
00443 }
00444 } else {
00445 return -1;
00446 }
00447
00448 cpl_errorstate status = cpl_errorstate_get();
00449
00450 cpl_table_new_column(aTable, "RA", CPL_TYPE_DOUBLE);
00451 cpl_table_new_column(aTable, "DEC", CPL_TYPE_DOUBLE);
00452
00453 cpl_size irow;
00454 for (irow = 0; irow < cpl_table_get_nrow(aTable); ++irow) {
00455 double x = cpl_table_get_double(aTable, "X", irow, NULL);
00456 double y = cpl_table_get_double(aTable, "Y", irow, NULL);
00457 double ra;
00458 double dec;
00459 muse_wcs_celestial_from_pixel_fast(wcs, x, y, &ra, &dec);
00460
00461 cpl_table_set_double(aTable, "RA", irow, ra);
00462 cpl_table_set_double(aTable, "DEC", irow, dec);
00463 }
00464 cpl_free(wcs);
00465
00466 if (!cpl_errorstate_is_equal(status)) {
00467 return -1;
00468 }
00469 return 0;
00470 }
00471
00472
00481
00482 static cpl_matrix *
00483 muse_cplmatrix_solve_least_square(const cpl_matrix *aMatrix1,
00484 const cpl_matrix *aMatrix2)
00485 {
00486 cpl_ensure(aMatrix1 && aMatrix2, CPL_ERROR_NULL_INPUT, NULL);
00487
00488 cpl_size nc = cpl_matrix_get_ncol(aMatrix1);
00489 cpl_size nr = cpl_matrix_get_nrow(aMatrix1);
00490 cpl_ensure(cpl_matrix_get_nrow(aMatrix2) == nr,
00491 CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
00492
00493 cpl_matrix *result = NULL;
00494 cpl_errorstate status = cpl_errorstate_get();
00495
00496 if (nr == nc) {
00497 result = cpl_matrix_solve(aMatrix1, aMatrix2);
00498 } else {
00499
00500 if (nr > nc) {
00501 result = cpl_matrix_solve_normal(aMatrix1, aMatrix2);
00502 } else {
00503 cpl_matrix *mt = cpl_matrix_transpose_create(aMatrix1);
00504 cpl_matrix *mgram = cpl_matrix_product_create(aMatrix1, mt);
00505 cpl_matrix *mw = cpl_matrix_solve(mgram, aMatrix2);
00506
00507 result = cpl_matrix_product_create(mt, mw);
00508 cpl_matrix_delete(mw);
00509 cpl_matrix_delete(mgram);
00510 cpl_matrix_delete(mt);
00511 }
00512 }
00513
00514 if (status != cpl_errorstate_get()) {
00515 cpl_matrix_delete(result);
00516 result = NULL;
00517 }
00518 return result;
00519 }
00520
00521
00531
00532 static cpl_error_code
00533 muse_cplmatrix_cosine(cpl_matrix *aMatrix)
00534 {
00535 cpl_ensure(aMatrix, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
00536
00537 cpl_size size = cpl_matrix_get_nrow(aMatrix) * cpl_matrix_get_ncol(aMatrix);
00538 double *mdata = cpl_matrix_get_data(aMatrix);
00539
00540 cpl_size i;
00541 for (i = 0; i < size; ++i) {
00542 mdata[i] = cos(mdata[i]);
00543 }
00544 return CPL_ERROR_NONE;
00545 }
00546
00547
00563
00564 static cpl_matrix *
00565 muse_align_compute_distances(const cpl_matrix *aDeltaRA, const cpl_matrix *aDeltaDEC,
00566 const cpl_matrix *aMeanDEC, const double *aOffsetRA,
00567 const double *aOffsetDEC)
00568 {
00569
00570 cpl_errorstate status = cpl_errorstate_get();
00571
00572
00573
00574
00575 cpl_matrix *dDEC;
00576 if (aOffsetDEC) {
00577 cpl_matrix *tmp = cpl_matrix_duplicate(aDeltaDEC);
00578 cpl_matrix_subtract_scalar(tmp, *aOffsetDEC);
00579
00580 dDEC = muse_cplmatrix_multiply_create(tmp, tmp);
00581 cpl_matrix_delete(tmp);
00582 } else {
00583 dDEC = muse_cplmatrix_multiply_create(aDeltaDEC, aDeltaDEC);
00584 }
00585
00586 cpl_matrix *dec = cpl_matrix_duplicate(aMeanDEC);
00587 muse_cplmatrix_cosine(dec);
00588
00589 cpl_matrix *dRA;
00590 if (aOffsetRA) {
00591 dRA = cpl_matrix_duplicate(aDeltaRA);
00592 cpl_matrix_subtract_scalar(dRA, *aOffsetRA);
00593 cpl_matrix_multiply(dRA, dec);
00594 } else {
00595 dRA = muse_cplmatrix_multiply_create(aDeltaRA, dec);
00596 }
00597 cpl_matrix_delete(dec);
00598
00599 cpl_matrix *distance = muse_cplmatrix_multiply_create(dRA, dRA);
00600 cpl_matrix_delete(dRA);
00601
00602 cpl_matrix_add(distance, dDEC);
00603 cpl_matrix_power(distance, 0.5);
00604 cpl_matrix_delete(dDEC);
00605
00606 if (status != cpl_errorstate_get()) {
00607 cpl_matrix_delete(distance);
00608 return NULL;
00609 }
00610
00611 return distance;
00612 }
00613
00614
00634
00635 static cpl_size
00636 muse_align_estimate_offsets(double *aOffsetRA, double *aOffsetDEC, double *aWeight,
00637 const cpl_matrix *aDeltaRA, const cpl_matrix *aDeltaDEC,
00638 const cpl_matrix *aMeanDEC, double aRadius, int aNbins)
00639 {
00640 cpl_ensure((aOffsetRA && aOffsetDEC) && aWeight, CPL_ERROR_NULL_INPUT, 0);
00641 cpl_ensure((aDeltaRA && aDeltaDEC) && aMeanDEC, CPL_ERROR_NULL_INPUT, 0);
00642 cpl_ensure(cpl_matrix_get_nrow(aDeltaRA) == cpl_matrix_get_nrow(aDeltaDEC),
00643 CPL_ERROR_INCOMPATIBLE_INPUT, 0);
00644 cpl_ensure(cpl_matrix_get_nrow(aDeltaRA) == cpl_matrix_get_nrow(aMeanDEC),
00645 CPL_ERROR_INCOMPATIBLE_INPUT, 0);
00646 cpl_ensure(aRadius > 0., CPL_ERROR_ILLEGAL_INPUT, 0);
00647 cpl_ensure(aNbins > 0, CPL_ERROR_ILLEGAL_INPUT, 0);
00648
00649 cpl_matrix *distance = muse_align_compute_distances(aDeltaRA, aDeltaDEC,
00650 aMeanDEC, NULL, NULL);
00651 if (cpl_matrix_get_min(distance) >= aRadius) {
00652 cpl_matrix_delete(distance);
00653 *aOffsetRA = 0.;
00654 *aOffsetDEC = 0.;
00655 *aWeight = 1.;
00656 return 0;
00657 }
00658
00659 cpl_array *selection = muse_cplmatrix_where(distance, aRadius,
00660 _muse_condition_less_than);
00661 cpl_matrix_delete(distance);
00662 cpl_matrix *dRA = muse_cplmatrix_extract_selected(aDeltaRA, selection);
00663 cpl_matrix *dDEC = muse_cplmatrix_extract_selected(aDeltaDEC, selection);
00664 cpl_array_delete(selection);
00665
00666
00667 double bstep = 2. * aRadius / (double)aNbins;
00668
00669
00670
00671
00672
00673
00674
00675
00676
00677
00678
00679
00680 cpl_matrix *bins = cpl_matrix_new(1, aNbins);
00681 cpl_size ibin;
00682 for (ibin = 0; ibin < aNbins; ++ibin) {
00683 cpl_matrix_set(bins, 0, ibin, -aRadius + ibin * bstep);
00684 }
00685 cpl_matrix *histogram = cpl_matrix_new(aNbins, aNbins);
00686 double *hdata = cpl_matrix_get_data(histogram);
00687 cpl_size npos = cpl_matrix_get_ncol(dRA);
00688 cpl_size nmatch = 0;
00689 cpl_size ipos;
00690 for (ipos = 0; ipos < npos; ++ipos) {
00691 cpl_size ix = (cpl_size)((cpl_matrix_get(dRA, 0, ipos) + aRadius) / bstep);
00692 cpl_size iy = (cpl_size)((cpl_matrix_get(dDEC, 0, ipos) + aRadius) / bstep);
00693 if (((ix >= 0) && (ix < aNbins)) && ((iy >= 0) && (iy < aNbins))) {
00694 hdata[iy * aNbins + ix] += 1.;
00695 ++nmatch;
00696 }
00697 }
00698 cpl_matrix_delete(dRA);
00699 cpl_matrix_delete(dDEC);
00700
00701 if (nmatch == 0) {
00702 cpl_matrix_delete(histogram);
00703 cpl_matrix_delete(bins);
00704 return 0;
00705 }
00706
00707 cpl_size row;
00708 cpl_size column;
00709 cpl_matrix_get_maxpos(histogram, &row, &column);
00710
00711 *aOffsetRA = cpl_matrix_get(bins, 0, column) + 0.5 * bstep;
00712 *aOffsetDEC = cpl_matrix_get(bins, 0, row) + 0.5 * bstep;
00713 *aWeight = cpl_matrix_get_max(histogram);
00714 cpl_matrix_delete(histogram);
00715 cpl_matrix_delete(bins);
00716
00717 return nmatch;
00718 }
00719
00720
00740
00741 static cpl_size
00742 muse_align_update_offsets(double *aOffsetRA, double *aOffsetDEC, double *aWeight,
00743 const cpl_matrix *aDeltaRA, const cpl_matrix *aDeltaDEC,
00744 const cpl_matrix *aMeanDEC, double aRadius)
00745 {
00746 cpl_ensure((aOffsetRA && aOffsetDEC) && aWeight, CPL_ERROR_NULL_INPUT,
00747 CPL_FALSE);
00748 cpl_ensure((aDeltaRA && aDeltaDEC) && aMeanDEC, CPL_ERROR_NULL_INPUT,
00749 CPL_FALSE);
00750 cpl_ensure(cpl_matrix_get_nrow(aDeltaRA) == cpl_matrix_get_nrow(aDeltaDEC),
00751 CPL_ERROR_INCOMPATIBLE_INPUT, CPL_FALSE);
00752 cpl_ensure(cpl_matrix_get_nrow(aDeltaRA) == cpl_matrix_get_nrow(aMeanDEC),
00753 CPL_ERROR_INCOMPATIBLE_INPUT, CPL_FALSE);
00754 cpl_ensure(aRadius > 0., CPL_ERROR_ILLEGAL_INPUT, CPL_FALSE);
00755
00756 cpl_matrix *distance = muse_align_compute_distances(aDeltaRA, aDeltaDEC,
00757 aMeanDEC, aOffsetRA, aOffsetDEC);
00758 if (cpl_matrix_get_min(distance) >= aRadius) {
00759 cpl_matrix_delete(distance);
00760 *aOffsetRA = 0.;
00761 *aOffsetDEC = 0.;
00762 *aWeight = 1.;
00763 return 0;
00764 }
00765
00766 cpl_array *selection = muse_cplmatrix_where(distance, aRadius,
00767 _muse_condition_less_than);
00768 cpl_matrix_delete(distance);
00769 cpl_matrix *dRA = muse_cplmatrix_extract_selected(aDeltaRA, selection);
00770 cpl_matrix *dDEC = muse_cplmatrix_extract_selected(aDeltaDEC, selection);
00771 cpl_array_delete(selection);
00772
00773 cpl_size nselected = cpl_matrix_get_ncol(dRA);
00774 double variance = 2. * aRadius * aRadius;
00775 if (nselected > 1) {
00776 double sdevRA = cpl_matrix_get_stdev(dRA);
00777 double sdevDEC = cpl_matrix_get_stdev(dDEC);
00778 double _variance = sdevRA * sdevRA + sdevDEC * sdevDEC;
00779 if (_variance > nselected / DBL_MAX) {
00780 variance = _variance;
00781 }
00782 }
00783
00784 *aOffsetRA = cpl_matrix_get_median(dRA);
00785 *aOffsetDEC = cpl_matrix_get_median(dDEC);
00786 *aWeight = nselected / variance;
00787
00788 cpl_matrix_delete(dRA);
00789 cpl_matrix_delete(dDEC);
00790
00791 return nselected;
00792 }
00793
00794
00809
00810 static cpl_table *
00811 muse_align_compute_field_offsets(muse_imagelist *aImagelist,
00812 cpl_table **aCataloglist,
00813 const cpl_array *aSearchradius,
00814 int aNbins, cpl_boolean aWeight,
00815 cpl_propertylist *aHeader)
00816 {
00817 cpl_ensure(aImagelist && aCataloglist, CPL_ERROR_NULL_INPUT, NULL);
00818 cpl_ensure(aSearchradius, CPL_ERROR_NULL_INPUT, NULL);
00819 cpl_ensure(aNbins > 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
00820
00821 cpl_size nfields = muse_imagelist_get_size(aImagelist);
00822 cpl_ensure(nfields >= 2, CPL_ERROR_DATA_NOT_FOUND, NULL);
00823
00824
00825
00826 cpl_size npairs = nfields * (nfields - 1) / 2;
00827
00828 muse_indexpair *combinations = cpl_calloc(npairs, sizeof *combinations);
00829 cpl_size ipair = 0;
00830 cpl_size ifield;
00831 for (ifield = 0; ifield < nfields; ++ifield) {
00832 cpl_size jfield;
00833 for (jfield = ifield + 1; jfield < nfields; ++jfield) {
00834 combinations[ipair].first = ifield;
00835 combinations[ipair].second = jfield;
00836 ++ipair;
00837 }
00838 }
00839
00840
00841 for (ifield = 0; ifield < nfields; ++ifield) {
00842 char *kw = cpl_sprintf(QC_ALIGN_NDETi, (int)ifield + 1);
00843 cpl_propertylist_append_int(aHeader, kw,
00844 cpl_table_get_nrow(aCataloglist[ifield]));
00845 cpl_free(kw);
00846 }
00847
00848
00849
00850
00851
00852 cpl_array **aoverlaps = cpl_calloc(nfields, sizeof(cpl_array *));
00853 for (ifield = 0; ifield < nfields; ++ifield) {
00854 aoverlaps[ifield] = cpl_array_new(nfields, CPL_TYPE_INT);
00855 }
00856
00857
00858 cpl_matrix *ra_offsets = NULL;
00859 cpl_matrix *dec_offsets = NULL;
00860 cpl_size *has_neighbors = cpl_malloc(nfields * sizeof *has_neighbors);
00861 cpl_size isearch;
00862 for (isearch = 0; isearch < cpl_array_get_size(aSearchradius); ++isearch) {
00863 double radius = cpl_array_get_double(aSearchradius, isearch, NULL);
00864
00865
00866
00867
00868
00869
00870
00871 for (ifield = 0; ifield < nfields; ++ifield) {
00872 cpl_array_fill_window_int(aoverlaps[ifield], 0, nfields, 0);
00873
00874
00875 #if 0
00876 cpl_array_set_invalid(aoverlaps[ifield], ifield);
00877 #else
00878 cpl_array_fill_window_invalid(aoverlaps[ifield], 0, nfields);
00879 #endif
00880 }
00881
00882
00883
00884 cpl_matrix *weights = cpl_matrix_new(npairs + 1, nfields);
00885 cpl_matrix *offsets_ra = cpl_matrix_new(npairs + 1, 1);
00886 cpl_matrix *offsets_dec = cpl_matrix_new(npairs + 1, 1);
00887
00888 cpl_matrix_fill_row(weights, 1., 0);
00889 cpl_size kpairs = 0;
00890
00891 memset(has_neighbors, 0, nfields * sizeof *has_neighbors);
00892
00893 for (ipair = 0; ipair < npairs; ++ipair) {
00894
00895 const cpl_table *catalog1 = aCataloglist[combinations[ipair].first];
00896 const cpl_table *catalog2 = aCataloglist[combinations[ipair].second];
00897
00898 cpl_size nstars1 = cpl_table_get_nrow(catalog1);
00899 cpl_size nstars2 = cpl_table_get_nrow(catalog2);
00900
00901 const double *ra1 = cpl_table_get_data_double_const(catalog1, "RA");
00902 const double *ra2 = cpl_table_get_data_double_const(catalog2, "RA");
00903 const double *dec1 = cpl_table_get_data_double_const(catalog1, "DEC");
00904 const double *dec2 = cpl_table_get_data_double_const(catalog2, "DEC");
00905
00906 cpl_matrix *delta_ra = cpl_matrix_new(nstars1, nstars2);
00907 cpl_matrix *delta_dec = cpl_matrix_new(nstars1, nstars2);
00908 cpl_matrix *dec_mean = cpl_matrix_new(nstars1, nstars2);
00909
00910 cpl_size irow;
00911 for (irow = 0; irow < nstars1; ++irow) {
00912 cpl_size icol;
00913 for (icol = 0; icol < nstars2; ++icol) {
00914 cpl_matrix_set(delta_ra, irow, icol, (ra1[irow] - ra2[icol]) * deg2as);
00915 cpl_matrix_set(delta_dec, irow, icol, (dec1[irow] - dec2[icol]) * deg2as);
00916 cpl_matrix_set(dec_mean, irow, icol,
00917 0.5 *(dec1[irow] + dec2[icol]) * CPL_MATH_RAD_DEG);
00918 }
00919 }
00920
00921 double offset_ra = 0.;
00922 double offset_dec = 0.;
00923 double weight = 1.;
00924 cpl_size noverlap = 0;
00925 if (isearch == 0) {
00926 noverlap = muse_align_estimate_offsets(&offset_ra, &offset_dec, &weight,
00927 delta_ra, delta_dec, dec_mean,
00928 radius, aNbins);
00929 } else {
00930 offset_ra = cpl_matrix_get(ra_offsets, combinations[ipair].first, 0) -
00931 cpl_matrix_get(ra_offsets, combinations[ipair].second, 0);
00932 offset_dec = cpl_matrix_get(dec_offsets, combinations[ipair].first, 0) -
00933 cpl_matrix_get(dec_offsets, combinations[ipair].second, 0);
00934
00935 noverlap = muse_align_update_offsets(&offset_ra, &offset_dec, &weight,
00936 delta_ra, delta_dec, dec_mean,
00937 radius);
00938 }
00939 cpl_matrix_delete(dec_mean);
00940 cpl_matrix_delete(delta_ra);
00941 cpl_matrix_delete(delta_dec);
00942
00943 if (noverlap > 0) {
00944 if (!aWeight) {
00945 weight = 1.;
00946 }
00947
00948 ++kpairs;
00949 cpl_matrix_set(offsets_ra, kpairs, 0, weight * offset_ra);
00950 cpl_matrix_set(offsets_dec, kpairs, 0, weight * offset_dec);
00951 cpl_matrix_set(weights, kpairs, combinations[ipair].first, weight);
00952 cpl_matrix_set(weights, kpairs, combinations[ipair].second, -weight);
00953 has_neighbors[combinations[ipair].first] += 1;
00954 has_neighbors[combinations[ipair].second] += 1;
00955
00956
00957 cpl_array_set_int(aoverlaps[combinations[ipair].first],
00958 combinations[ipair].second, noverlap);
00959 cpl_array_set_int(aoverlaps[combinations[ipair].second],
00960 combinations[ipair].first, noverlap);
00961 }
00962 }
00963
00964
00965 cpl_matrix_delete(ra_offsets);
00966 cpl_matrix_delete(dec_offsets);
00967
00968
00969
00970 if (!kpairs) {
00971 cpl_matrix_delete(offsets_dec);
00972 cpl_matrix_delete(offsets_ra);
00973 cpl_matrix_delete(weights);
00974
00975 cpl_msg_warning(__func__, "No overlapping fields of view were found "
00976 "for search radius %.4f!", radius);
00977
00978 ra_offsets = cpl_matrix_new(nfields, 1);
00979 dec_offsets = cpl_matrix_new(nfields, 1);
00980
00981 } else {
00982
00983 cpl_matrix *_weights = cpl_matrix_wrap(kpairs + 1, nfields,
00984 cpl_matrix_get_data(weights));
00985 cpl_matrix *_offsets_ra = cpl_matrix_wrap(kpairs + 1, 1,
00986 cpl_matrix_get_data(offsets_ra));
00987 cpl_matrix *_offsets_dec = cpl_matrix_wrap(kpairs + 1, 1,
00988 cpl_matrix_get_data(offsets_dec));
00989
00990 ra_offsets = muse_cplmatrix_solve_least_square(_weights, _offsets_ra);
00991 dec_offsets = muse_cplmatrix_solve_least_square(_weights, _offsets_dec);
00992
00993 cpl_matrix_unwrap(_offsets_dec);
00994 cpl_matrix_unwrap(_offsets_ra);
00995 cpl_matrix_unwrap(_weights);
00996
00997 cpl_matrix_delete(offsets_dec);
00998 cpl_matrix_delete(offsets_ra);
00999 cpl_matrix_delete(weights);
01000 }
01001 }
01002 cpl_free(combinations);
01003
01004
01005
01006
01007
01008
01009 cpl_table *offsets = muse_cpltable_new(muse_offset_list_def,
01010 muse_imagelist_get_size(aImagelist));
01011 cpl_size minmatch = LLONG_MAX,
01012 nomatch = 0;
01013 for (ifield = 0; ifield < nfields; ++ifield) {
01014 muse_image *fov = muse_imagelist_get(aImagelist, ifield);
01015 const char *timestamp = muse_pfits_get_dateobs(fov->header);
01016 double mjd = muse_pfits_get_mjdobs(fov->header);
01017 cpl_table_set_string(offsets, MUSE_OFFSETS_DATEOBS, ifield, timestamp);
01018 cpl_table_set_double(offsets, MUSE_OFFSETS_MJDOBS, ifield, mjd);
01019
01020 if (has_neighbors[ifield]) {
01021 double ra_offset = cpl_matrix_get(ra_offsets, ifield, 0);
01022 double dec_offset = cpl_matrix_get(dec_offsets, ifield, 0);
01023 cpl_table_set_double(offsets, MUSE_OFFSETS_DRA, ifield, ra_offset / deg2as);
01024 cpl_table_set_double(offsets, MUSE_OFFSETS_DDEC, ifield, dec_offset / deg2as);
01025
01026 int nmedian = cpl_array_get_median(aoverlaps[ifield]);
01027 char *kw = cpl_sprintf(QC_ALIGN_NMATCHi, (int)ifield + 1);
01028 cpl_propertylist_append_int(aHeader, kw, nmedian);
01029 cpl_free(kw);
01030 if (nmedian < minmatch) {
01031 minmatch = nmedian;
01032 }
01033 } else {
01034 cpl_table_set_double(offsets, MUSE_OFFSETS_DRA, ifield, 0.);
01035 cpl_table_set_double(offsets, MUSE_OFFSETS_DDEC, ifield, 0.);
01036 char *kw = cpl_sprintf(QC_ALIGN_NMATCHi, (int)ifield + 1);
01037 cpl_propertylist_append_int(aHeader, kw, 0);
01038 cpl_free(kw);
01039 nomatch++;
01040 }
01041 }
01042 cpl_propertylist_append_int(aHeader, QC_ALIGN_NMATCH_MIN, minmatch);
01043 cpl_propertylist_append_int(aHeader, QC_ALIGN_NOMATCH, nomatch);
01044 muse_vfree((void **)aoverlaps, nfields, (muse_free_func)cpl_array_delete);
01045 cpl_free(has_neighbors);
01046 cpl_matrix_delete(dec_offsets);
01047 cpl_matrix_delete(ra_offsets);
01048
01049
01050
01051
01052 double ra_shift = cpl_table_get(offsets, MUSE_OFFSETS_DRA, 0, NULL);
01053 double dec_shift = cpl_table_get(offsets, MUSE_OFFSETS_DDEC, 0, NULL);
01054 cpl_table_subtract_scalar(offsets, MUSE_OFFSETS_DRA, ra_shift);
01055 cpl_table_subtract_scalar(offsets, MUSE_OFFSETS_DDEC, dec_shift);
01056 cpl_table_set_column_invalid(offsets, MUSE_OFFSETS_FSCALE, 0, nfields);
01057
01058
01059 double qcmin = cpl_table_get_column_min(offsets, MUSE_OFFSETS_DRA);
01060 double qcmax = cpl_table_get_column_max(offsets, MUSE_OFFSETS_DRA);
01061 double qcmean = cpl_table_get_column_mean(offsets, MUSE_OFFSETS_DRA);
01062 double qcsdev = cpl_table_get_column_stdev(offsets, MUSE_OFFSETS_DRA);
01063
01064 cpl_propertylist_append_double(aHeader, QC_ALIGN_DRA_MIN, qcmin);
01065 cpl_propertylist_append_double(aHeader, QC_ALIGN_DRA_MAX, qcmax);
01066 cpl_propertylist_append_double(aHeader, QC_ALIGN_DRA_MEAN, qcmean);
01067 cpl_propertylist_append_double(aHeader, QC_ALIGN_DRA_STDEV, qcsdev);
01068
01069 qcmin = cpl_table_get_column_min(offsets, MUSE_OFFSETS_DDEC);
01070 qcmax = cpl_table_get_column_max(offsets, MUSE_OFFSETS_DDEC);
01071 qcmean = cpl_table_get_column_mean(offsets, MUSE_OFFSETS_DDEC);
01072 qcsdev = cpl_table_get_column_stdev(offsets, MUSE_OFFSETS_DDEC);
01073
01074 cpl_propertylist_append_double(aHeader, QC_ALIGN_DDEC_MIN, qcmin);
01075 cpl_propertylist_append_double(aHeader, QC_ALIGN_DDEC_MAX, qcmax);
01076 cpl_propertylist_append_double(aHeader, QC_ALIGN_DDEC_MEAN, qcmean);
01077 cpl_propertylist_append_double(aHeader, QC_ALIGN_DDEC_STDEV, qcsdev);
01078
01079 return offsets;
01080 }
01081
01082
01091
01092 int
01093 muse_exp_align_compute(muse_processing *aProcessing,
01094 muse_exp_align_params_t *aParams)
01095 {
01096
01097 cpl_size nexposures = cpl_frameset_count_tags(aProcessing->inframes,
01098 MUSE_TAG_IMAGE_FOV);
01099 if (nexposures < 2) {
01100 cpl_msg_error(__func__, "This recipe requires at least 2 input FOV "
01101 "images!");
01102 return -1;
01103 }
01104
01105
01106 cpl_msg_debug(__func__, "Validating source detection parameters");
01107 if (!muse_align_check_detection_params(aParams)) {
01108 cpl_error_set_where(__func__);
01109 return -1;
01110 }
01111
01112
01113 cpl_msg_info(__func__, "Loading FOV images");
01114
01115 muse_imagelist *fovimages = muse_processing_fov_load_all(aProcessing);
01116 if (!fovimages) {
01117 cpl_msg_error(__func__, "Could not load FOV input images!");
01118 return -1;
01119 }
01120 nexposures = (cpl_size)muse_imagelist_get_size(fovimages);
01121
01122 cpl_size iexposure;
01123 for (iexposure = 0; iexposure < nexposures; ++iexposure) {
01124 const muse_image *fov = muse_imagelist_get(fovimages, iexposure);
01125 cpl_errorstate es = cpl_errorstate_get();
01126 double mjdobs = muse_pfits_get_mjdobs(fov->header);
01127
01128 if (!cpl_errorstate_is_equal(es) &&
01129 cpl_error_get_code() == CPL_ERROR_DATA_NOT_FOUND) {
01130 cpl_msg_error(__func__, "Missing \"MJD-OBS\" in FOV image %"
01131 CPL_SIZE_FORMAT, iexposure + 1);
01132 muse_imagelist_delete(fovimages);
01133 return -1;
01134 } else {
01135 cpl_size jexposure;
01136 for (jexposure = iexposure + 1; jexposure < nexposures; ++jexposure) {
01137 const muse_image *_fov = muse_imagelist_get(fovimages, jexposure);
01138 double _mjdobs = muse_pfits_get_mjdobs(_fov->header);
01139 if (_mjdobs == mjdobs) {
01140 cpl_msg_warning(__func__, "Found FOV images with non-unique "
01141 "\"MJD-OBS\" value (image %" CPL_SIZE_FORMAT
01142 " and %" CPL_SIZE_FORMAT ")!", iexposure + 1,
01143 jexposure + 1);
01144 }
01145 }
01146 }
01147 }
01148
01149 cpl_msg_info(__func__, "Computing pointing corrections for %" CPL_SIZE_FORMAT
01150 " FOV images", nexposures);
01151
01152
01153
01154 for (iexposure = 0; iexposure < nexposures; ++iexposure) {
01155 muse_image *image = muse_imagelist_get(fovimages, iexposure);
01156 muse_utils_replace_nan(image, 0.);
01157 }
01158
01159
01160
01161 int nobjects[2] = {aParams->srcmin, aParams->srcmax};
01162
01163 double roundness[2] = {aParams->roundmin, aParams->roundmax};
01164 double sharpness[2] = {aParams->sharpmin, aParams->sharpmax};
01165
01166 cpl_table **srclists = cpl_calloc(nexposures, sizeof *srclists);
01167
01168 for (iexposure = 0; iexposure < nexposures; ++iexposure) {
01169 cpl_msg_debug(__func__, "Detecting sources on image %" CPL_SIZE_FORMAT
01170 " of %" CPL_SIZE_FORMAT, iexposure + 1, nexposures);
01171
01172 muse_image *fov = muse_imagelist_get(fovimages, iexposure);
01173 double threshold = aParams->threshold;
01174 cpl_boolean iterate = CPL_TRUE;
01175 int iteration = 0;
01176 cpl_size ndetections = 0;
01177 cpl_table *detections;
01178 cpl_image *image = cpl_image_cast(fov->data, CPL_TYPE_DOUBLE);
01179
01180
01181
01182
01183
01184
01185
01186 cpl_size nx = cpl_image_get_size_x(image);
01187 cpl_size ny = cpl_image_get_size_y(image);
01188
01189 if ((nx < 6) || (ny < 6)) {
01190 cpl_image_delete(image);
01191 muse_vfree((void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
01192 muse_imagelist_delete(fovimages);
01193 cpl_msg_error(__func__, "Image %" CPL_SIZE_FORMAT " of %"
01194 CPL_SIZE_FORMAT "is too small!", iexposure + 1, nexposures);
01195 return -1;
01196 }
01197
01198 if ((nx % 2 != 0) || (ny % 2 != 0)) {
01199 cpl_msg_debug(__func__, "Trimming image %" CPL_SIZE_FORMAT " to an "
01200 "even number of pixels along the axes.", iexposure + 1);
01201 if (nx % 2 != 0) {
01202 --nx;
01203 }
01204 if (ny % 2 != 0) {
01205 --ny;
01206 }
01207 cpl_image *_image = cpl_image_extract(image, 1, 1, nx, ny);
01208 if (!_image) {
01209 cpl_image_delete(image);
01210 muse_vfree((void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
01211 muse_imagelist_delete(fovimages);
01212 cpl_msg_error(__func__, "Trimming image %" CPL_SIZE_FORMAT " failed",
01213 iexposure + 1);
01214 return -1;
01215 }
01216 cpl_image_delete(image);
01217 image = _image;
01218 }
01219
01220 while (iterate) {
01221 cpl_errorstate status = cpl_errorstate_get();
01222
01223 detections = muse_find_stars(image, threshold, aParams->fwhm,
01224 roundness, sharpness);
01225 if (!detections) {
01226 if (cpl_error_get_code() == CPL_ERROR_DATA_NOT_FOUND) {
01227 cpl_errorstate_set(status);
01228 ndetections = 0;
01229 } else {
01230 cpl_image_delete(image);
01231 muse_vfree((void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
01232 muse_imagelist_delete(fovimages);
01233 cpl_msg_error(__func__, "Source detection failed on image %"
01234 CPL_SIZE_FORMAT " of %" CPL_SIZE_FORMAT, iexposure + 1,
01235 nexposures);
01236 return -1;
01237 }
01238 } else {
01239 ndetections = cpl_table_get_nrow(detections);
01240 }
01241
01242 if ((ndetections > 0) &&
01243 ((ndetections >= nobjects[0]) && (ndetections <= nobjects[1]))) {
01244 iterate = CPL_FALSE;
01245 } else {
01246 cpl_table_delete(detections);
01247 detections = NULL;
01248 if (++iteration < aParams->iterations) {
01249 if (ndetections > nobjects[1]) {
01250 threshold += aParams->step;
01251 } else {
01252 threshold -= aParams->step;
01253 if (threshold <= 0.) {
01254 iterate = CPL_FALSE;
01255 }
01256 }
01257 } else {
01258 iterate = CPL_FALSE;
01259 }
01260 }
01261 }
01262 cpl_image_delete(image);
01263
01264 if (!detections) {
01265 muse_vfree((void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
01266 muse_imagelist_delete(fovimages);
01267 if (ndetections == 0) {
01268 cpl_msg_error(__func__, "No point-sources were detected in image %"
01269 CPL_SIZE_FORMAT, iexposure + 1);
01270 } else if (ndetections < nobjects[0]) {
01271 cpl_msg_error(__func__, "Too few (%" CPL_SIZE_FORMAT ") point-sources "
01272 "were detected in image %" CPL_SIZE_FORMAT, ndetections,
01273 iexposure + 1);
01274 } else if (ndetections > nobjects[1]) {
01275 cpl_msg_error(__func__, "Too many (%" CPL_SIZE_FORMAT ") point-sources "
01276 "were detected in image %" CPL_SIZE_FORMAT, ndetections,
01277 iexposure + 1);
01278 } else {
01279 cpl_msg_error(__func__, "Invalid number (%" CPL_SIZE_FORMAT ") of "
01280 "point-sources were detected in image %" CPL_SIZE_FORMAT,
01281 ndetections, iexposure + 1);
01282 }
01283 return -1;
01284 }
01285
01286 cpl_msg_info(__func__, "%" CPL_SIZE_FORMAT " point-sources detected in "
01287 "image %" CPL_SIZE_FORMAT, ndetections, iexposure + 1);
01288
01289
01290
01291 if (muse_align_celestial_from_pixel(detections, fov->header)) {
01292 muse_vfree((void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
01293 muse_imagelist_delete(fovimages);
01294 cpl_msg_error(__func__, "Computing WCS coordinates from pixel positions "
01295 "failed for image %" CPL_SIZE_FORMAT " of %"
01296 CPL_SIZE_FORMAT, iexposure + 1, nexposures);
01297 return -1;
01298 }
01299 srclists[iexposure] = detections;
01300 }
01301
01302
01303
01304 cpl_array *radii = muse_align_parse_valuelist(aParams->rsearch);
01305 if (!radii) {
01306 muse_vfree((void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
01307 muse_imagelist_delete(fovimages);
01308 cpl_msg_error(__func__, "Cannot compute field offsets! No valid search "
01309 "radius was given!");
01310 return -1;
01311 }
01312
01313 cpl_msg_info(__func__, "Calculating field offsets...");
01314 cpl_propertylist *header = cpl_propertylist_new();
01315 cpl_table *offsets = muse_align_compute_field_offsets(fovimages, srclists, radii,
01316 aParams->nbins, aParams->weight,
01317 header);
01318 cpl_array_delete(radii);
01319 muse_vfree((void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
01320 muse_imagelist_delete(fovimages);
01321 if (!offsets) {
01322 cpl_msg_error(__func__, "Could not compute FOV offsets for %"
01323 CPL_SIZE_FORMAT " images", nexposures);
01324 cpl_propertylist_delete(header);
01325 return -1;
01326 }
01327
01328
01329 cpl_error_code rc = muse_processing_save_table(aProcessing, -1, offsets,
01330 header, MUSE_TAG_OFFSET_LIST,
01331 MUSE_TABLE_TYPE_CPL);
01332 cpl_propertylist_delete(header);
01333 cpl_table_delete(offsets);
01334
01335 if (rc != CPL_ERROR_NONE) {
01336 return -1;
01337 }
01338
01339 return 0;
01340 }