MUSE Pipeline Reference Manual  2.1.1
muse_exp_align.c
1 /* -*- Mode: C; tab-width: 2; indent-tabs-mode: nil; c-basic-offset: 2 -*- */
2 /* vim:set sw=2 sts=2 et cin: */
3 /*
4  * This file is part of the MUSE Instrument Pipeline
5  * Copyright (C) 2015 European Southern Observatory
6  *
7  * This program is free software; you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation; either version 2 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this program; if not, write to the Free Software
19  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
20  */
21 
22 #ifdef HAVE_CONFIG_H
23 #include <config.h>
24 #endif
25 
26 /*---------------------------------------------------------------------------*
27  * Includes *
28  *---------------------------------------------------------------------------*/
29 #include <errno.h>
30 #include <stdlib.h>
31 #include <string.h>
32 #include <strings.h>
33 #include <math.h>
34 #include <cpl.h>
35 #include <muse.h>
36 #include "muse_exp_align_z.h"
37 
38 
39 typedef struct muse_indexpair muse_indexpair;
40 
41 struct muse_indexpair {
42  cpl_size first;
43  cpl_size second;
44 };
45 
46 typedef void (*muse_free_func)(void *);
47 
48 
49 /* Conversion factor degree to arcsecond */
50 static const double deg2as = 3600.;
51 
52 
53 /*---------------------------------------------------------------------------*
54  * Functions code *
55  *---------------------------------------------------------------------------*/
56 
57 static void
58 muse_vfree(void **array, cpl_size size, muse_free_func deallocator)
59 {
60  if (array) {
61  cpl_size idx;
62  for (idx = 0; idx < size; ++idx) {
63  if (deallocator) {
64  deallocator(array[idx]);
65  }
66  }
67  cpl_free(array);
68  }
69  return;
70 }
71 
72 static int
73 _muse_condition_less_than(double aValue1, double aValue2)
74 {
75  return (aValue1 < aValue2) ? TRUE : FALSE;
76 }
77 
78 /*----------------------------------------------------------------------------*/
92 /*----------------------------------------------------------------------------*/
93 static cpl_error_code
94 muse_utils_replace_nan(muse_image *aImage, float aValue)
95 {
96  cpl_ensure(aImage, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
97  cpl_ensure(aImage->data, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_ILLEGAL_INPUT);
98 
99  cpl_size npixel = cpl_image_get_size_x(aImage->data) *
100  cpl_image_get_size_y(aImage->data);
101 
102  float *data = cpl_image_get_data_float(aImage->data);
103  if (cpl_error_get_code() == CPL_ERROR_TYPE_MISMATCH) {
104  cpl_error_set_message(__func__, CPL_ERROR_TYPE_MISMATCH,
105  "MUSE image data buffer has invalid type!");
106  return CPL_ERROR_TYPE_MISMATCH;
107  }
108 
109  cpl_size ipixel;
110  for (ipixel = 0; ipixel < npixel; ++ipixel) {
111  if (isnan(data[ipixel])) {
112  data[ipixel] = aValue;
113  }
114  }
115 
116  if (aImage->stat) {
117  npixel = cpl_image_get_size_x(aImage->stat) *
118  cpl_image_get_size_y(aImage->stat);
119 
120  data = cpl_image_get_data_float(aImage->stat);
121  if (cpl_error_get_code() == CPL_ERROR_TYPE_MISMATCH) {
122  cpl_error_set_message(__func__, CPL_ERROR_TYPE_MISMATCH,
123  "MUSE image data buffer has invalid type!");
124  return CPL_ERROR_TYPE_MISMATCH;
125  }
126 
127  for (ipixel = 0; ipixel < npixel; ++ipixel) {
128  if (isnan(data[ipixel])) {
129  data[ipixel] = aValue;
130  }
131  }
132  }
133 
134  return CPL_ERROR_NONE;
135 }
136 
137 /* XXX: Moved to muse_utils, kept here disabled until final place is found */
138 #if 0
139 /*----------------------------------------------------------------------------*/
155 /*----------------------------------------------------------------------------*/
156 static muse_image *
157 muse_fov_load(const char *aFilename)
158 {
159  muse_image *image = muse_image_new();
160 
161  image->header = cpl_propertylist_load(aFilename, 0);
162  if (!image->header) {
163  cpl_error_set_message(__func__, cpl_error_get_code(), "Loading primary FITS "
164  "header of \"%s\" did not succeed", aFilename);
165  muse_image_delete(image);
166  return NULL;
167  }
168 
169  /* assuming an IMAGE_FOV input, start by searching for the DATA extension */
170  int extension = cpl_fits_find_extension(aFilename, EXTNAME_DATA);
171  cpl_propertylist *hdata = cpl_propertylist_load(aFilename, extension);
172  while (muse_pfits_get_naxis(hdata, 0) != 2) {
173  /* Not an image, this means we were given a cube, possibly *
174  * with image extensions; search for the first of those. */
175  cpl_msg_debug(__func__, "Skipping extension %d [%s]", extension,
176  muse_pfits_get_extname(hdata));
177  cpl_propertylist_delete(hdata);
178  extension++;
179  hdata = cpl_propertylist_load(aFilename, extension);
180  } /* while */
181  cpl_msg_debug(__func__, "Taking extension %d [%s]", extension,
182  muse_pfits_get_extname(hdata));
183  /* keep the extension name, to check if we loaded *
184  * from an IMAGE_FOV or from a cube extension */
185  char *extname = cpl_strdup(muse_pfits_get_extname(hdata));
186  image->data = cpl_image_load(aFilename, CPL_TYPE_FLOAT, 0, extension);
187  if (!image->data) {
188  cpl_error_set_message(__func__, MUSE_ERROR_READ_DATA, "Could not load "
189  "extension %s from \"%s\"", extname, aFilename);
190  muse_image_delete(image);
191  cpl_free(extname);
192  return NULL;
193  }
194 
195  if (cpl_propertylist_has(hdata, "BUNIT")) {
196  cpl_propertylist_append_string(image->header, "BUNIT",
197  muse_pfits_get_bunit(hdata));
198  cpl_propertylist_set_comment(image->header, "BUNIT",
199  cpl_propertylist_get_comment(hdata, "BUNIT"));
200  } else {
201  cpl_msg_warning(__func__, "No BUNIT given in extension %d [%s] of \"%s\"!",
202  extension, extname, aFilename);
203  }
204 
205  if (!cpl_propertylist_has(hdata, "CUNIT1") ||
206  !cpl_propertylist_has(hdata, "CUNIT2")) {
207  cpl_msg_warning(__func__, "No WCS found in extension %d [%s] of \"%s\"!",
208  extension, extname, aFilename);
209  }
210 
211  /* Merge WCS keywords and ESO hierarchical keywords from the data *
212  * extension into the image header */
213  cpl_propertylist_erase_regexp(hdata, "(^ESO |" MUSE_WCS_KEYS ")", 1);
214  cpl_propertylist_append(image->header, hdata);
215  cpl_propertylist_delete(hdata);
216 
217  /* Load STAT extension if it is found and contains data, ignore it *
218  * otherwise */
219  if (extname && !strncmp(extname, EXTNAME_DATA, strlen(EXTNAME_DATA)+1)) {
220  /* if this is an IMAGE_FOV, then a STAT extension might be there */
221  extension = cpl_fits_find_extension(aFilename, EXTNAME_STAT);
222  } else {
223  /* otherwise we might have an extension <filtername>_STAT */
224  char *statname = cpl_sprintf("%s_STAT", extname);
225  extension = cpl_fits_find_extension(aFilename, statname);
226  cpl_free(statname);
227  }
228  if (extension > 0) {
229  cpl_errorstate status = cpl_errorstate_get();
230  image->stat = cpl_image_load(aFilename, CPL_TYPE_INT, 0, extension);
231  if (!cpl_errorstate_is_equal(status)) {
232  if (cpl_error_get_code() == CPL_ERROR_DATA_NOT_FOUND) {
233  cpl_errorstate_set(status);
234  cpl_msg_debug(__func__, "Ignoring empty extension %s in \"%s\"",
235  EXTNAME_STAT, aFilename);
236  } else {
237  cpl_error_set_message(__func__, MUSE_ERROR_READ_STAT, "Could not load "
238  "extension %s from \"%s\"", EXTNAME_STAT, aFilename);
239  muse_image_delete(image);
240  cpl_free(extname);
241  return NULL;
242  }
243  }
244  }
245 
246  /* If present load the DQ extension and encode flagged pixels in the *
247  * image data (and possibly the statistics as NaNs. */
248  if (extname && !strncmp(extname, EXTNAME_DATA, strlen(EXTNAME_DATA)+1)) {
249  /* if this is an IMAGE_FOV, then a DQ extension might be there */
250  extension = cpl_fits_find_extension(aFilename, EXTNAME_DQ);
251  } else {
252  /* otherwise we might have an extension <filtername>_DQ */
253  char *dqname = cpl_sprintf("%s_DQ", extname);
254  extension = cpl_fits_find_extension(aFilename, dqname);
255  cpl_free(dqname);
256  }
257  if (extension > 0) {
258  cpl_errorstate status = cpl_errorstate_get();
259  image->dq = cpl_image_load(aFilename, CPL_TYPE_INT, 0, extension);
260  if (!cpl_errorstate_is_equal(status)) {
261  cpl_errorstate_set(status);
262  cpl_error_set_message(__func__, MUSE_ERROR_READ_DQ, "Could not load "
263  "extension %s from \"%s\"", EXTNAME_DQ, aFilename);
264  muse_image_delete(image);
265  cpl_free(extname);
266  return NULL;
267  }
268  muse_image_dq_to_nan(image);
269  }
270 
271  cpl_free(extname);
272  return image;
273 }
274 #endif
275 
276 /* XXX: Maybe there is a better name for this function */
277 /*----------------------------------------------------------------------------*/
283 /*----------------------------------------------------------------------------*/
284 static muse_imagelist *
285 muse_processing_fov_load_all(muse_processing *aProcessing)
286 {
287  cpl_ensure(aProcessing, CPL_ERROR_NULL_INPUT, NULL);
288  cpl_size nframes = cpl_frameset_get_size(aProcessing->inframes);
289  cpl_ensure(nframes, CPL_ERROR_DATA_NOT_FOUND, NULL);
290 
291  muse_imagelist *fovimages = muse_imagelist_new();
292 
293  cpl_size iexposure = 0;
294  cpl_size iframe;
295  for (iframe = 0; iframe < nframes; ++iframe) {
296  const cpl_frame *frame = cpl_frameset_get_position(aProcessing->inframes,
297  iframe);
298  const char *tag = cpl_frame_get_tag(frame);
299 
300  if (muse_processing_check_intags(aProcessing, tag, strlen(tag))) {
301  const char *filename = cpl_frame_get_filename(frame);
302 
303  cpl_msg_debug(__func__, "Loading FOV image '%s' as exposure %"
304  CPL_SIZE_FORMAT, filename, iexposure + 1);
305 
306  muse_image *image = muse_fov_load(filename);
307  if (!image) {
308  cpl_msg_error(__func__, "Could not load FOV image '%s'", filename);
309  muse_imagelist_delete(fovimages);
310  return NULL;
311  }
312 
313  muse_imagelist_set(fovimages, image, iexposure);
314  muse_processing_append_used(aProcessing, (cpl_frame *)frame,
315  CPL_FRAME_GROUP_RAW, 1);
316  ++iexposure;
317  }
318  }
319  return fovimages;
320 }
321 
322 /*----------------------------------------------------------------------------*/
329 /*----------------------------------------------------------------------------*/
330 static cpl_boolean
331 muse_align_check_detection_params(muse_exp_align_params_t *aParams)
332 {
333  if (aParams->threshold <= 0.) {
334  cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
335  "Detection threshold must be greater than zero!");
336  return CPL_FALSE;
337  }
338  if (aParams->fwhm < 0.5) {
339  cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
340  "FWHM must be greater than 0.5 pixels!");
341  return CPL_FALSE;
342  }
343  if (aParams->srcmax <= aParams->srcmin) {
344  cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
345  "Maximum number of sources must be grater than the "
346  "minimum number of sources!");
347  return CPL_FALSE;
348  }
349  if (aParams->roundmax < aParams->roundmin) {
350  cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
351  "Upper limit of the point-source roundness must be "
352  "greater than the lower limit!");
353  return CPL_FALSE;
354  }
355  if (aParams->sharpmin <= 0.) {
356  cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
357  "Low limit of the point-source sharpness must be "
358  "greater than zero!");
359  return CPL_FALSE;
360  }
361  if (aParams->sharpmax < aParams->sharpmin) {
362  cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
363  "Upper limit of the point-source sharpness must be "
364  "greater than the lower limit!");
365  return CPL_FALSE;
366  }
367  return CPL_TRUE;
368 }
369 
370 /*----------------------------------------------------------------------------*/
378 /*----------------------------------------------------------------------------*/
379 static cpl_array *
380 muse_align_parse_valuelist(const char *aValuelist)
381 {
382  cpl_ensure(aValuelist, CPL_ERROR_NULL_INPUT, NULL);
383  if (strlen(aValuelist) == 0) {
384  cpl_error_set_message(__func__, CPL_ERROR_DATA_NOT_FOUND,
385  "List of values is empty!");
386  return NULL;
387  }
388 
389  cpl_array *strings = muse_cplarray_new_from_delimited_string(aValuelist, ",");
390  cpl_size nvalues = cpl_array_get_size(strings);
391 
392  cpl_array *values = cpl_array_new(nvalues, CPL_TYPE_DOUBLE);
393  cpl_size ivalue;
394  for (ivalue = 0; ivalue < nvalues; ++ivalue) {
395  const char *svalue = cpl_array_get_string(strings, ivalue);
396  if ((strncasecmp(svalue, "inf", 3) == 0) ||
397  (strncasecmp(svalue, "nan", 3) == 0)) {
398  cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
399  "Illegal value \"%s\" encountered!", svalue);
400  cpl_array_delete(values);
401  cpl_array_delete(strings);
402  return NULL;
403  }
404  else {
405  char *last = NULL;
406  double value = strtod(svalue, &last);
407  if( (errno == ERANGE) || ((value == 0.) && (last == svalue))) {
408  cpl_error_set_message(__func__, CPL_ERROR_TYPE_MISMATCH,
409  "Could not convert \"%s\" to a decimal number!", svalue);
410  cpl_array_delete(values);
411  cpl_array_delete(strings);
412  return NULL;
413  }
414  cpl_array_set_double(values, ivalue, value);
415  }
416  }
417  cpl_array_delete(strings);
418  return values;
419 }
420 
421 /*----------------------------------------------------------------------------*/
430 /*----------------------------------------------------------------------------*/
431 static int
432 muse_align_celestial_from_pixel(cpl_table *aTable,
433  cpl_propertylist *aWCSHeader)
434 {
435  muse_wcs *wcs = muse_wcs_new(aWCSHeader);
436 
437  const char *unit1 = muse_pfits_get_cunit(aWCSHeader, 1);
438  const char *unit2 = muse_pfits_get_cunit(aWCSHeader, 2);
439 
440  if (unit1 && unit2) {
441  if ((strncmp(unit1, unit2, 4) == 0) && (strncmp(unit1, "deg", 3) == 0)) {
442  wcs->iscelsph = CPL_TRUE;
443  }
444  } else {
445  return -1;
446  }
447 
448  cpl_errorstate status = cpl_errorstate_get();
449 
450  cpl_table_new_column(aTable, "RA", CPL_TYPE_DOUBLE);
451  cpl_table_new_column(aTable, "DEC", CPL_TYPE_DOUBLE);
452 
453  cpl_size irow;
454  for (irow = 0; irow < cpl_table_get_nrow(aTable); ++irow) {
455  double x = cpl_table_get_double(aTable, "X", irow, NULL);
456  double y = cpl_table_get_double(aTable, "Y", irow, NULL);
457  double ra;
458  double dec;
459  muse_wcs_celestial_from_pixel_fast(wcs, x, y, &ra, &dec);
460 
461  cpl_table_set_double(aTable, "RA", irow, ra);
462  cpl_table_set_double(aTable, "DEC", irow, dec);
463  }
464  cpl_free(wcs);
465 
466  if (!cpl_errorstate_is_equal(status)) {
467  return -1;
468  }
469  return 0;
470 }
471 
472 /* XXX: Keep the conditional until the SVD code has been fully verified. *
473  * Once this is done, completely replace the old code with the SVD *
474  * version. */
475 #define USE_SVD_SOLVER
476 #ifdef USE_SVD_SOLVER
477 
478 #if CPL_VERSION_CODE < 459008
479 #define cpl_matrix_solve_svd muse_cplmatrix_solve_svd
480 
481 /*
482  * The following implementation of one-sided Jacobi orthogonalization
483  * to compute the SVD of a MxN matrix, with M >= N is based on the
484  * implementation found in the GNU Scientific Library (GSL).
485  */
486 
487 /*----------------------------------------------------------------------------*/
521 /*----------------------------------------------------------------------------*/
522 static cpl_error_code
523 muse_cplmatrix_decomp_sv_jacobi(cpl_matrix *U, cpl_vector *S, cpl_matrix *V)
524 {
525  cpl_ensure_code((U != NULL) && (V != NULL) && (S != NULL),
526  CPL_ERROR_NULL_INPUT);
527  cpl_ensure_code(cpl_matrix_get_nrow(U) >= cpl_matrix_get_ncol(U),
528  CPL_ERROR_ILLEGAL_INPUT);
529  cpl_ensure_code(cpl_matrix_get_nrow(V) == cpl_matrix_get_ncol(V),
530  CPL_ERROR_ILLEGAL_INPUT);
531  cpl_ensure_code(cpl_matrix_get_nrow(V) == cpl_matrix_get_ncol(U),
532  CPL_ERROR_INCOMPATIBLE_INPUT);
533  cpl_ensure_code(cpl_vector_get_size(S) == cpl_matrix_get_ncol(U),
534  CPL_ERROR_INCOMPATIBLE_INPUT);
535 
536  const cpl_size m = cpl_matrix_get_nrow(U);
537  const cpl_size n = cpl_matrix_get_ncol(U);
538 
539  /* Initialize the rotation counter and the sweep counter. *
540  * Use a minimum number of 12 sweeps. */
541  cpl_size count = 1;
542  cpl_size sweep = 0;
543  cpl_size sweepmax = CPL_MAX(5 * n, 12);
544 
545  double tolerance = 10. * m * DBL_EPSILON;
546 
547  /* Fill V with the identity matrix. */
548  cpl_matrix_fill(V, 0.);
549  cpl_matrix_fill_diagonal(V, 1., 0);
550 
551  /* Store the column error estimates in S, for use during the *
552  * orthogonalization. */
553 
554  cpl_size i;
555  cpl_size j;
556  cpl_vector *cv1 = cpl_vector_new(m);
557  double *const _cv1 = cpl_vector_get_data(cv1);
558 
559  for (j = 0; j < n; ++j) {
560  for (i = 0; i < m; ++i) {
561  _cv1[i] = cpl_matrix_get(U, i, j);
562  }
563 
564  double s = sqrt(cpl_vector_product(cv1, cv1));
565  cpl_vector_set(S, j, DBL_EPSILON * s);
566  }
567 
568  /* Orthogonalization of U by plane rotations. */
569  cpl_vector *cv2 = cpl_vector_new(m);
570  double *const _cv2 = cpl_vector_get_data(cv2);
571 
572  while ((count > 0) && (sweep <= sweepmax)) {
573 
574  /* Initialize the rotation counter. */
575  count = n * (n - 1) / 2;
576 
577  for (j = 0; j < n - 1; ++j) {
578  cpl_size k;
579  for (k = j + 1; k < n; ++k) {
580  for (i = 0; i < m; ++i) {
581  _cv1[i] = cpl_matrix_get(U, i, j);
582  }
583  for (i = 0; i < m; ++i) {
584  _cv2[i] = cpl_matrix_get(U, i, k);
585  }
586 
587  /* Compute the (j, k) submatrix elements of Ut*U. The *
588  * non-diagonal element c is scaled by 2. */
589  double a = sqrt(cpl_vector_product(cv1, cv1));
590  double b = sqrt(cpl_vector_product(cv2, cv2));
591  double c = 2. * cpl_vector_product(cv1, cv2);
592 
593  /* Test whether columns j, k are orthogonal, or dominant errors. */
594  double abserr_a = cpl_vector_get(S, j);
595  double abserr_b = cpl_vector_get(S, k);
596  cpl_boolean orthogonal = (fabs(c) <= tolerance * a * b);
597  cpl_boolean sorted = (a >= b);
598  cpl_boolean noisy_a = (a < abserr_a);
599  cpl_boolean noisy_b = (b < abserr_b);
600 
601  if (sorted && (orthogonal || noisy_a || noisy_b)) {
602  --count;
603  continue;
604  }
605 
606  /* Compute the Jacobi rotation which diagonalizes the (j, k) *
607  * submatrix. */
608  double q = a * a - b * b;
609  double v = sqrt(c * c + q * q);
610  double cs; /* cos(theta) */
611  double sn; /* sin(theta) */
612 
613  if (v == 0. || !sorted) {
614  cs = 0.;
615  sn = 1.;
616  } else {
617  cs = sqrt((v + q) / (2. * v));
618  sn = c / (2. * v * cs);
619  }
620 
621  /* Update the columns j and k of U and V */
622  for (i = 0; i < m; ++i) {
623  const double u_ij = cpl_matrix_get(U, i, j);
624  const double u_ik = cpl_matrix_get(U, i, k);
625  cpl_matrix_set(U, i, j, u_ij * cs + u_ik * sn);
626  cpl_matrix_set(U, i, k, -u_ij * sn + u_ik * cs);
627  }
628 
629  for (i = 0; i < n; ++i) {
630  const double v_ij = cpl_matrix_get(V, i, j);
631  const double v_ik = cpl_matrix_get(V, i, k);
632  cpl_matrix_set(V, i, j, v_ij * cs + v_ik * sn);
633  cpl_matrix_set(V, i, k, -v_ij * sn + v_ik * cs);
634  }
635 
636  cpl_vector_set(S, j, fabs(cs) * abserr_a +
637  fabs(sn) * abserr_b);
638  cpl_vector_set(S, k, fabs(sn) * abserr_a +
639  fabs(cs) * abserr_b);
640  }
641  }
642  ++sweep;
643 
644  }
645  cpl_vector_delete(cv2);
646 
647  /* Compute singular values. */
648  double prev_norm = -1.;
649  for (j = 0; j < n; ++j) {
650  for (i = 0; i < m; ++i) {
651  _cv1[i] = cpl_matrix_get(U, i, j);
652  }
653 
654  double norm = sqrt(cpl_vector_product(cv1, cv1));
655 
656  /* Determine if the singular value is zero, according to the *
657  * criteria used in the main loop above (i.e. comparison with *
658  * norm of previous column). */
659  if ((norm == 0.) || (prev_norm == 0.) ||
660  ((j > 0) && (norm <= tolerance * prev_norm))) {
661  cpl_vector_set(S, j, 0.); /* singular */
662  cpl_matrix_fill_column(U, 0., j); /* annihilate column */
663  prev_norm = 0.;
664  } else {
665  cpl_vector_set(S, j, norm); /* non-singular */
666 
667  /* Normalize the column vector of U */
668  for (i = 0; i < m; ++i) {
669  cpl_matrix_set(U, i, j, _cv1[i] / norm);
670  }
671  prev_norm = norm;
672  }
673  }
674  cpl_vector_delete(cv1);
675 
676  if (count > 0) {
677  return CPL_ERROR_CONTINUE;
678  }
679  return CPL_ERROR_NONE;
680 }
681 
682 /*----------------------------------------------------------------------------*/
691 /*----------------------------------------------------------------------------*/
692 static cpl_matrix *
693 muse_cplmatrix_solve_svd(const cpl_matrix *coeff, const cpl_matrix *rhs)
694 {
695  cpl_ensure((coeff != NULL) && (rhs != NULL), CPL_ERROR_NULL_INPUT, NULL);
696  cpl_ensure(cpl_matrix_get_nrow(coeff) == cpl_matrix_get_nrow(rhs),
697  CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
698  cpl_ensure(cpl_matrix_get_ncol(rhs) == 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
699 
700  /* Compute the SVD factorization of the coefficient matrix */
701  cpl_matrix *U = cpl_matrix_duplicate(coeff);
702  cpl_size nc = cpl_matrix_get_ncol(U);
703 
704  cpl_matrix *V = cpl_matrix_new(nc, nc);
705  cpl_vector *S = cpl_vector_new(nc);
706 
707  cpl_error_code status = muse_cplmatrix_decomp_sv_jacobi(U, S, V);
708  if (status != CPL_ERROR_NONE) {
709  cpl_vector_delete(S);
710  cpl_matrix_delete(V);
711  cpl_matrix_delete(U);
712  cpl_error_set_where(cpl_func);
713  return NULL;
714  }
715 
716  /* Solve the linear system */
717  cpl_matrix *Ut = cpl_matrix_transpose_create(U);
718  cpl_matrix *w = cpl_matrix_product_create(Ut, rhs);
719  cpl_matrix_delete(Ut);
720 
721  const double *_s = cpl_vector_get_data_const(S);
722  double *_w = cpl_matrix_get_data(w);
723  cpl_size ic;
724 
725  for (ic = 0; ic < nc; ++ic) {
726  double s = _s[ic];
727  if (s != 0.) {
728  s = 1. / s;
729  }
730  _w[ic] *= s;
731  }
732 
733  cpl_matrix *solution = cpl_matrix_product_create(V, w);
734  cpl_matrix_delete(w);
735  cpl_vector_delete(S);
736  cpl_matrix_delete(V);
737  cpl_matrix_delete(U);
738  return solution;
739 }
740 
741 #endif
742 
743 #endif
744 
745 /*----------------------------------------------------------------------------*/
754 /*----------------------------------------------------------------------------*/
755 static cpl_matrix *
756 muse_cplmatrix_solve_least_square(const cpl_matrix *aMatrix1,
757  const cpl_matrix *aMatrix2)
758 {
759  cpl_ensure(aMatrix1 && aMatrix2, CPL_ERROR_NULL_INPUT, NULL);
760 
761  cpl_size nc = cpl_matrix_get_ncol(aMatrix1);
762  cpl_size nr = cpl_matrix_get_nrow(aMatrix1);
763  cpl_ensure(cpl_matrix_get_nrow(aMatrix2) == nr,
764  CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
765 
766  cpl_matrix *result = NULL;
767  cpl_errorstate status = cpl_errorstate_get();
768 
769 #ifndef USE_SVD_SOLVER
770  if (nr == nc) {
771  result = cpl_matrix_solve(aMatrix1, aMatrix2);
772  } else {
773  /* Handle over- and under-determined systems */
774  if (nr > nc) {
775  result = cpl_matrix_solve_normal(aMatrix1, aMatrix2);
776  } else {
777  cpl_matrix *mt = cpl_matrix_transpose_create(aMatrix1);
778  cpl_matrix *mgram = cpl_matrix_product_create(aMatrix1, mt);
779  cpl_matrix *mw = cpl_matrix_solve(mgram, aMatrix2);
780 
781  result = cpl_matrix_product_create(mt, mw);
782  cpl_matrix_delete(mw);
783  cpl_matrix_delete(mgram);
784  cpl_matrix_delete(mt);
785  }
786  }
787 #else
788  if (nr >= nc) {
789  result = cpl_matrix_solve_svd(aMatrix1, aMatrix2);
790  } else {
791  cpl_matrix *mt = cpl_matrix_transpose_create(aMatrix1);
792  cpl_matrix *mgram = cpl_matrix_product_create(aMatrix1, mt);
793  cpl_matrix *mw = cpl_matrix_solve(mgram, aMatrix2);
794 
795  result = cpl_matrix_product_create(mt, mw);
796  cpl_matrix_delete(mw);
797  cpl_matrix_delete(mgram);
798  cpl_matrix_delete(mt);
799  }
800 #endif
801 
802  if (status != cpl_errorstate_get()) {
803  cpl_matrix_delete(result);
804  result = NULL;
805  }
806  return result;
807 }
808 
809 /*----------------------------------------------------------------------------*/
819 /*----------------------------------------------------------------------------*/
820 static cpl_error_code
821 muse_cplmatrix_cosine(cpl_matrix *aMatrix)
822 {
823  cpl_ensure(aMatrix, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
824 
825  cpl_size size = cpl_matrix_get_nrow(aMatrix) * cpl_matrix_get_ncol(aMatrix);
826  double *mdata = cpl_matrix_get_data(aMatrix);
827 
828  cpl_size i;
829  for (i = 0; i < size; ++i) {
830  mdata[i] = cos(mdata[i]);
831  }
832  return CPL_ERROR_NONE;
833 }
834 
835 /*----------------------------------------------------------------------------*/
851 /*----------------------------------------------------------------------------*/
852 static cpl_matrix *
853 muse_align_compute_distances(const cpl_matrix *aDeltaRA, const cpl_matrix *aDeltaDEC,
854  const cpl_matrix *aMeanDEC, const double *aOffsetRA,
855  const double *aOffsetDEC)
856 {
857 
858  cpl_errorstate status = cpl_errorstate_get();
859 
860  /* XXX: Cannot use cpl_matrix_power() here, since it chokes on negative *
861  * elements, however it should not if the exponent is integer! */
862 
863  cpl_matrix *dDEC;
864  if (aOffsetDEC) {
865  cpl_matrix *tmp = cpl_matrix_duplicate(aDeltaDEC);
866  cpl_matrix_subtract_scalar(tmp, *aOffsetDEC);
867 
868  dDEC = muse_cplmatrix_multiply_create(tmp, tmp);
869  cpl_matrix_delete(tmp);
870  } else {
871  dDEC = muse_cplmatrix_multiply_create(aDeltaDEC, aDeltaDEC);
872  }
873 
874  cpl_matrix *dec = cpl_matrix_duplicate(aMeanDEC);
875  muse_cplmatrix_cosine(dec);
876 
877  cpl_matrix *dRA;
878  if (aOffsetRA) {
879  dRA = cpl_matrix_duplicate(aDeltaRA);
880  cpl_matrix_subtract_scalar(dRA, *aOffsetRA);
881  cpl_matrix_multiply(dRA, dec);
882  } else {
883  dRA = muse_cplmatrix_multiply_create(aDeltaRA, dec);
884  }
885  cpl_matrix_delete(dec);
886 
887  cpl_matrix *distance = muse_cplmatrix_multiply_create(dRA, dRA);
888  cpl_matrix_delete(dRA);
889 
890  cpl_matrix_add(distance, dDEC);
891  cpl_matrix_power(distance, 0.5);
892  cpl_matrix_delete(dDEC);
893 
894  if (status != cpl_errorstate_get()) {
895  cpl_matrix_delete(distance);
896  return NULL;
897  }
898 
899  return distance;
900 }
901 
902 /*----------------------------------------------------------------------------*/
922 /*----------------------------------------------------------------------------*/
923 static cpl_size
924 muse_align_estimate_offsets(double *aOffsetRA, double *aOffsetDEC, double *aWeight,
925  const cpl_matrix *aDeltaRA, const cpl_matrix *aDeltaDEC,
926  const cpl_matrix *aMeanDEC, double aRadius, int aNbins)
927 {
928  cpl_ensure((aOffsetRA && aOffsetDEC) && aWeight, CPL_ERROR_NULL_INPUT, 0);
929  cpl_ensure((aDeltaRA && aDeltaDEC) && aMeanDEC, CPL_ERROR_NULL_INPUT, 0);
930  cpl_ensure(cpl_matrix_get_nrow(aDeltaRA) == cpl_matrix_get_nrow(aDeltaDEC),
931  CPL_ERROR_INCOMPATIBLE_INPUT, 0);
932  cpl_ensure(cpl_matrix_get_nrow(aDeltaRA) == cpl_matrix_get_nrow(aMeanDEC),
933  CPL_ERROR_INCOMPATIBLE_INPUT, 0);
934  cpl_ensure(aRadius > 0., CPL_ERROR_ILLEGAL_INPUT, 0);
935  cpl_ensure(aNbins > 0, CPL_ERROR_ILLEGAL_INPUT, 0);
936 
937  cpl_matrix *distance = muse_align_compute_distances(aDeltaRA, aDeltaDEC,
938  aMeanDEC, NULL, NULL);
939  if (cpl_matrix_get_min(distance) >= aRadius) {
940  cpl_matrix_delete(distance);
941  *aOffsetRA = 0.;
942  *aOffsetDEC = 0.;
943  *aWeight = 1.;
944  return 0;
945  }
946 
947  cpl_array *selection = muse_cplmatrix_where(distance, aRadius,
948  _muse_condition_less_than);
949  cpl_matrix_delete(distance);
950  cpl_matrix *dRA = muse_cplmatrix_extract_selected(aDeltaRA, selection);
951  cpl_matrix *dDEC = muse_cplmatrix_extract_selected(aDeltaDEC, selection);
952  cpl_array_delete(selection);
953 
954  /* Search box from -aRadius to aRadius for both axes. */
955  double bstep = 2. * aRadius / (double)aNbins;
956 
957  /* Note: Despite of the previous selection of positions which have *
958  * distances smaller than the search radius the bin index of the right *
959  * ascension may become invalid and must be restricted to valid range. *
960  * The reason is that in the following the difference in right ascension *
961  * is directly compared to the search radius, while in the selection of *
962  * positions this difference is corrected by the cosine of the mean *
963  * declination (cf. muse_align_compute_distances()), and thus differences *
964  * in right ascension which are larger than the search radius may still *
965  * be accepted by the distance criterion! *
966  * To play it save also the index along the y-axis is validated. */
967 
968  cpl_matrix *bins = cpl_matrix_new(1, aNbins);
969  cpl_size ibin;
970  for (ibin = 0; ibin < aNbins; ++ibin) {
971  cpl_matrix_set(bins, 0, ibin, -aRadius + ibin * bstep);
972  }
973  cpl_matrix *histogram = cpl_matrix_new(aNbins, aNbins);
974  double *hdata = cpl_matrix_get_data(histogram);
975  cpl_size npos = cpl_matrix_get_ncol(dRA);
976  cpl_size nmatch = 0;
977  cpl_size ipos;
978  for (ipos = 0; ipos < npos; ++ipos) {
979  cpl_size ix = (cpl_size)((cpl_matrix_get(dRA, 0, ipos) + aRadius) / bstep);
980  cpl_size iy = (cpl_size)((cpl_matrix_get(dDEC, 0, ipos) + aRadius) / bstep);
981  if (((ix >= 0) && (ix < aNbins)) && ((iy >= 0) && (iy < aNbins))) {
982  hdata[iy * aNbins + ix] += 1.;
983  ++nmatch;
984  }
985  }
986  cpl_matrix_delete(dRA);
987  cpl_matrix_delete(dDEC);
988 
989  if (nmatch == 0) {
990  cpl_matrix_delete(histogram);
991  cpl_matrix_delete(bins);
992  return 0;
993  }
994 
995  cpl_size row;
996  cpl_size column;
997  cpl_matrix_get_maxpos(histogram, &row, &column);
998 
999  *aOffsetRA = cpl_matrix_get(bins, 0, column) + 0.5 * bstep;
1000  *aOffsetDEC = cpl_matrix_get(bins, 0, row) + 0.5 * bstep;
1001  *aWeight = cpl_matrix_get_max(histogram);
1002  cpl_matrix_delete(histogram);
1003  cpl_matrix_delete(bins);
1004 
1005  return nmatch;
1006 }
1007 
1008 /*----------------------------------------------------------------------------*/
1028 /*----------------------------------------------------------------------------*/
1029 static cpl_size
1030 muse_align_update_offsets(double *aOffsetRA, double *aOffsetDEC, double *aWeight,
1031  const cpl_matrix *aDeltaRA, const cpl_matrix *aDeltaDEC,
1032  const cpl_matrix *aMeanDEC, double aRadius)
1033 {
1034  cpl_ensure((aOffsetRA && aOffsetDEC) && aWeight, CPL_ERROR_NULL_INPUT,
1035  CPL_FALSE);
1036  cpl_ensure((aDeltaRA && aDeltaDEC) && aMeanDEC, CPL_ERROR_NULL_INPUT,
1037  CPL_FALSE);
1038  cpl_ensure(cpl_matrix_get_nrow(aDeltaRA) == cpl_matrix_get_nrow(aDeltaDEC),
1039  CPL_ERROR_INCOMPATIBLE_INPUT, CPL_FALSE);
1040  cpl_ensure(cpl_matrix_get_nrow(aDeltaRA) == cpl_matrix_get_nrow(aMeanDEC),
1041  CPL_ERROR_INCOMPATIBLE_INPUT, CPL_FALSE);
1042  cpl_ensure(aRadius > 0., CPL_ERROR_ILLEGAL_INPUT, CPL_FALSE);
1043 
1044  cpl_matrix *distance = muse_align_compute_distances(aDeltaRA, aDeltaDEC,
1045  aMeanDEC, aOffsetRA, aOffsetDEC);
1046  if (cpl_matrix_get_min(distance) >= aRadius) {
1047  cpl_matrix_delete(distance);
1048  *aOffsetRA = 0.;
1049  *aOffsetDEC = 0.;
1050  *aWeight = 1.;
1051  return 0;
1052  }
1053 
1054  cpl_array *selection = muse_cplmatrix_where(distance, aRadius,
1055  _muse_condition_less_than);
1056  cpl_matrix_delete(distance);
1057  cpl_matrix *dRA = muse_cplmatrix_extract_selected(aDeltaRA, selection);
1058  cpl_matrix *dDEC = muse_cplmatrix_extract_selected(aDeltaDEC, selection);
1059  cpl_array_delete(selection);
1060 
1061  cpl_size nselected = cpl_matrix_get_ncol(dRA);
1062  double variance = 2. * aRadius * aRadius;
1063  if (nselected > 1) {
1064  double sdevRA = cpl_matrix_get_stdev(dRA);
1065  double sdevDEC = cpl_matrix_get_stdev(dDEC);
1066  double _variance = sdevRA * sdevRA + sdevDEC * sdevDEC;
1067  if (_variance > nselected / DBL_MAX) {
1068  variance = _variance;
1069  }
1070  }
1071 
1072  *aOffsetRA = cpl_matrix_get_median(dRA);
1073  *aOffsetDEC = cpl_matrix_get_median(dDEC);
1074  *aWeight = nselected / variance;
1075 
1076  cpl_matrix_delete(dRA);
1077  cpl_matrix_delete(dDEC);
1078 
1079  return nselected;
1080 }
1081 
1082 /*----------------------------------------------------------------------------*/
1097 /*----------------------------------------------------------------------------*/
1098 static cpl_table *
1099 muse_align_compute_field_offsets(muse_imagelist *aImagelist,
1100  cpl_table **aCataloglist,
1101  const cpl_array *aSearchradius,
1102  int aNbins, cpl_boolean aWeight,
1103  cpl_propertylist *aHeader)
1104 {
1105  cpl_ensure(aImagelist && aCataloglist, CPL_ERROR_NULL_INPUT, NULL);
1106  cpl_ensure(aSearchradius, CPL_ERROR_NULL_INPUT, NULL);
1107  cpl_ensure(aNbins > 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
1108 
1109  cpl_size nfields = muse_imagelist_get_size(aImagelist);
1110  cpl_ensure(nfields >= 2, CPL_ERROR_DATA_NOT_FOUND, NULL);
1111 
1112  /* Create a list of all possible pairings of the input fields *
1113  * number of combinations: nfields choose 2 */
1114  cpl_size npairs = nfields * (nfields - 1) / 2;
1115 
1116  muse_indexpair *combinations = cpl_calloc(npairs, sizeof *combinations);
1117  cpl_size ipair = 0;
1118  cpl_size ifield;
1119  for (ifield = 0; ifield < nfields; ++ifield) {
1120  cpl_size jfield;
1121  for (jfield = ifield + 1; jfield < nfields; ++jfield) {
1122  combinations[ipair].first = ifield;
1123  combinations[ipair].second = jfield;
1124  ++ipair;
1125  }
1126  }
1127 
1128  /* record number of detections in each image as QC parameter */
1129  for (ifield = 0; ifield < nfields; ++ifield) {
1130  char *kw = cpl_sprintf(QC_ALIGN_NDETi, (int)ifield + 1);
1131  cpl_propertylist_append_int(aHeader, kw,
1132  cpl_table_get_nrow(aCataloglist[ifield]));
1133  cpl_free(kw);
1134  }
1135 
1136  /* prepare the data structure to record the number of matching *
1137  * objects found in each of the other fields to, eventually, *
1138  * compute the median number of matching detections for each *
1139  * field. */
1140  cpl_array **aoverlaps = cpl_calloc(nfields, sizeof(cpl_array *));
1141  for (ifield = 0; ifield < nfields; ++ifield) {
1142  aoverlaps[ifield] = cpl_array_new(nfields, CPL_TYPE_INT);
1143  }
1144 
1145  /* search for correlations */
1146  cpl_matrix *ra_offsets = NULL;
1147  cpl_matrix *dec_offsets = NULL;
1148  cpl_size *has_neighbors = cpl_malloc(nfields * sizeof *has_neighbors);
1149  cpl_size isearch;
1150  for (isearch = 0; isearch < cpl_array_get_size(aSearchradius); ++isearch) {
1151  double radius = cpl_array_get_double(aSearchradius, isearch, NULL);
1152 
1153  /* Clear the number of matching detections. Note that this resets all *
1154  * elements to valid. To exclude the elements which should not be taken *
1155  * into account in the final median computation they are invalidated. *
1156  * This must at least be done for the "diagonal" elements of this *
1157  * matrix like structure which would refer to the number of matching *
1158  * detections a field has with itself, since these are never computed. */
1159  for (ifield = 0; ifield < nfields; ++ifield) {
1160  cpl_array_fill_window_int(aoverlaps[ifield], 0, nfields, 0);
1161  /* XXX: To include combinations with zero matches in the median *
1162  * statistic, use the alternative invalidation scheme */
1163 #if 0
1164  cpl_array_set_invalid(aoverlaps[ifield], ifield);
1165 #else
1166  cpl_array_fill_window_invalid(aoverlaps[ifield], 0, nfields);
1167 #endif
1168  }
1169 
1170  /* The linear system consists of one row for each pair of fields, and *
1171  * the closure relation */
1172  cpl_matrix *weights = cpl_matrix_new(npairs + 1, nfields);
1173  cpl_matrix *offsets_ra = cpl_matrix_new(npairs + 1, 1);
1174  cpl_matrix *offsets_dec = cpl_matrix_new(npairs + 1, 1);
1175 
1176  cpl_matrix_fill_row(weights, 1., 0);
1177  cpl_size kpairs = 0;
1178 
1179  memset(has_neighbors, 0, nfields * sizeof *has_neighbors);
1180 
1181  for (ipair = 0; ipair < npairs; ++ipair) {
1182  /* Compute the distance between any 2 stars in the 2 fields */
1183  const cpl_table *catalog1 = aCataloglist[combinations[ipair].first];
1184  const cpl_table *catalog2 = aCataloglist[combinations[ipair].second];
1185 
1186  cpl_size nstars1 = cpl_table_get_nrow(catalog1);
1187  cpl_size nstars2 = cpl_table_get_nrow(catalog2);
1188 
1189  const double *ra1 = cpl_table_get_data_double_const(catalog1, "RA");
1190  const double *ra2 = cpl_table_get_data_double_const(catalog2, "RA");
1191  const double *dec1 = cpl_table_get_data_double_const(catalog1, "DEC");
1192  const double *dec2 = cpl_table_get_data_double_const(catalog2, "DEC");
1193 
1194  cpl_matrix *delta_ra = cpl_matrix_new(nstars1, nstars2);
1195  cpl_matrix *delta_dec = cpl_matrix_new(nstars1, nstars2);
1196  cpl_matrix *dec_mean = cpl_matrix_new(nstars1, nstars2);
1197 
1198  cpl_size irow;
1199  for (irow = 0; irow < nstars1; ++irow) {
1200  cpl_size icol;
1201  for (icol = 0; icol < nstars2; ++icol) {
1202  cpl_matrix_set(delta_ra, irow, icol, (ra1[irow] - ra2[icol]) * deg2as);
1203  cpl_matrix_set(delta_dec, irow, icol, (dec1[irow] - dec2[icol]) * deg2as);
1204  cpl_matrix_set(dec_mean, irow, icol,
1205  0.5 *(dec1[irow] + dec2[icol]) * CPL_MATH_RAD_DEG);
1206  }
1207  }
1208 
1209  double offset_ra = 0.;
1210  double offset_dec = 0.;
1211  double weight = 1.;
1212  cpl_size noverlap = 0;
1213  if (isearch == 0) {
1214  noverlap = muse_align_estimate_offsets(&offset_ra, &offset_dec, &weight,
1215  delta_ra, delta_dec, dec_mean,
1216  radius, aNbins);
1217  } else {
1218  offset_ra = cpl_matrix_get(ra_offsets, combinations[ipair].first, 0) -
1219  cpl_matrix_get(ra_offsets, combinations[ipair].second, 0);
1220  offset_dec = cpl_matrix_get(dec_offsets, combinations[ipair].first, 0) -
1221  cpl_matrix_get(dec_offsets, combinations[ipair].second, 0);
1222 
1223  noverlap = muse_align_update_offsets(&offset_ra, &offset_dec, &weight,
1224  delta_ra, delta_dec, dec_mean,
1225  radius);
1226  }
1227  cpl_matrix_delete(dec_mean);
1228  cpl_matrix_delete(delta_ra);
1229  cpl_matrix_delete(delta_dec);
1230 
1231  if (noverlap > 0) {
1232  if (!aWeight) {
1233  weight = 1.;
1234  }
1235 
1236  ++kpairs;
1237  cpl_matrix_set(offsets_ra, kpairs, 0, weight * offset_ra);
1238  cpl_matrix_set(offsets_dec, kpairs, 0, weight * offset_dec);
1239  cpl_matrix_set(weights, kpairs, combinations[ipair].first, weight);
1240  cpl_matrix_set(weights, kpairs, combinations[ipair].second, -weight);
1241  has_neighbors[combinations[ipair].first] += 1;
1242  has_neighbors[combinations[ipair].second] += 1;
1243 
1244  /* record the number of overlaps for later median computation */
1245  cpl_array_set_int(aoverlaps[combinations[ipair].first],
1246  combinations[ipair].second, noverlap);
1247  cpl_array_set_int(aoverlaps[combinations[ipair].second],
1248  combinations[ipair].first, noverlap);
1249  } /* if noverlap positive */
1250  } /* for ipair */
1251 
1252  /* Replace the solution from the previous iteration */
1253  cpl_matrix_delete(ra_offsets);
1254  cpl_matrix_delete(dec_offsets);
1255 
1256  /* If there is no field with overlapping neighbors there is no *
1257  * correction to the field coordinates, i.e. all offsets are 0. */
1258  if (!kpairs) {
1259  cpl_matrix_delete(offsets_dec);
1260  cpl_matrix_delete(offsets_ra);
1261  cpl_matrix_delete(weights);
1262 
1263  cpl_msg_warning(__func__, "No overlapping fields of view were found "
1264  "for search radius %.4f!", radius);
1265 
1266  ra_offsets = cpl_matrix_new(nfields, 1);
1267  dec_offsets = cpl_matrix_new(nfields, 1);
1268 
1269  } else {
1270  /* Use only the valid part of the matrices when solving the linear system */
1271  cpl_matrix *_weights = cpl_matrix_wrap(kpairs + 1, nfields,
1272  cpl_matrix_get_data(weights));
1273  cpl_matrix *_offsets_ra = cpl_matrix_wrap(kpairs + 1, 1,
1274  cpl_matrix_get_data(offsets_ra));
1275  cpl_matrix *_offsets_dec = cpl_matrix_wrap(kpairs + 1, 1,
1276  cpl_matrix_get_data(offsets_dec));
1277 
1278  ra_offsets = muse_cplmatrix_solve_least_square(_weights, _offsets_ra);
1279  dec_offsets = muse_cplmatrix_solve_least_square(_weights, _offsets_dec);
1280 
1281  cpl_matrix_unwrap(_offsets_dec);
1282  cpl_matrix_unwrap(_offsets_ra);
1283  cpl_matrix_unwrap(_weights);
1284 
1285  cpl_matrix_delete(offsets_dec);
1286  cpl_matrix_delete(offsets_ra);
1287  cpl_matrix_delete(weights);
1288  }
1289  } /* for isearch */
1290  cpl_free(combinations);
1291 
1292  /* Build the results table from the matrices, and add observation start *
1293  * from the image header as observation identifier. For isolated fields, *
1294  * i.e.fields which do not overlap with any other, the pointing *
1295  * corrections are forced to 0. For fields which have an offset could be *
1296  * computed it is converted to degrees. */
1297  cpl_table *offsets = muse_cpltable_new(muse_offset_list_def,
1298  muse_imagelist_get_size(aImagelist));
1299  cpl_size minmatch = LLONG_MAX,
1300  nomatch = 0;
1301  for (ifield = 0; ifield < nfields; ++ifield) {
1302  muse_image *fov = muse_imagelist_get(aImagelist, ifield);
1303  const char *timestamp = muse_pfits_get_dateobs(fov->header);
1304  double mjd = muse_pfits_get_mjdobs(fov->header);
1305  cpl_table_set_string(offsets, MUSE_OFFSETS_DATEOBS, ifield, timestamp);
1306  cpl_table_set_double(offsets, MUSE_OFFSETS_MJDOBS, ifield, mjd);
1307 
1308  if (has_neighbors[ifield]) {
1309  double ra_offset = cpl_matrix_get(ra_offsets, ifield, 0);
1310  double dec_offset = cpl_matrix_get(dec_offsets, ifield, 0);
1311  cpl_table_set_double(offsets, MUSE_OFFSETS_DRA, ifield, ra_offset / deg2as);
1312  cpl_table_set_double(offsets, MUSE_OFFSETS_DDEC, ifield, dec_offset / deg2as);
1313  /* record QC parameter */
1314  int nmedian = cpl_array_get_median(aoverlaps[ifield]);
1315  char *kw = cpl_sprintf(QC_ALIGN_NMATCHi, (int)ifield + 1);
1316  cpl_propertylist_append_int(aHeader, kw, nmedian);
1317  cpl_free(kw);
1318  if (nmedian < minmatch) {
1319  minmatch = nmedian;
1320  }
1321  } else {
1322  cpl_table_set_double(offsets, MUSE_OFFSETS_DRA, ifield, 0.);
1323  cpl_table_set_double(offsets, MUSE_OFFSETS_DDEC, ifield, 0.);
1324  char *kw = cpl_sprintf(QC_ALIGN_NMATCHi, (int)ifield + 1);
1325  cpl_propertylist_append_int(aHeader, kw, 0);
1326  cpl_free(kw);
1327  nomatch++;
1328  }
1329  } /* for ifield */
1330  cpl_propertylist_append_int(aHeader, QC_ALIGN_NMATCH_MIN, minmatch);
1331  cpl_propertylist_append_int(aHeader, QC_ALIGN_NOMATCH, nomatch);
1332  muse_vfree((void **)aoverlaps, nfields, (muse_free_func)cpl_array_delete);
1333  cpl_free(has_neighbors);
1334  cpl_matrix_delete(dec_offsets);
1335  cpl_matrix_delete(ra_offsets);
1336 
1337  /* Force the offset of the first field of view to be zero and *
1338  * invalidate the flux scaling column, so that by no scaling *
1339  * occurs during the exposure combination. */
1340  double ra_shift = cpl_table_get(offsets, MUSE_OFFSETS_DRA, 0, NULL);
1341  double dec_shift = cpl_table_get(offsets, MUSE_OFFSETS_DDEC, 0, NULL);
1342  cpl_table_subtract_scalar(offsets, MUSE_OFFSETS_DRA, ra_shift);
1343  cpl_table_subtract_scalar(offsets, MUSE_OFFSETS_DDEC, dec_shift);
1344  cpl_table_set_column_invalid(offsets, MUSE_OFFSETS_FSCALE, 0, nfields);
1345 
1346  /* Extra QC parameters: offset statistics, in arcsec */
1347  double qcmin = cpl_table_get_column_min(offsets, MUSE_OFFSETS_DRA) * 3600.;
1348  double qcmax = cpl_table_get_column_max(offsets, MUSE_OFFSETS_DRA) * 3600.;
1349  double qcmean = cpl_table_get_column_mean(offsets, MUSE_OFFSETS_DRA) * 3600.;
1350  double qcsdev = cpl_table_get_column_stdev(offsets, MUSE_OFFSETS_DRA) * 3600.;
1351 
1352  cpl_propertylist_append_double(aHeader, QC_ALIGN_DRA_MIN, qcmin);
1353  cpl_propertylist_append_double(aHeader, QC_ALIGN_DRA_MAX, qcmax);
1354  cpl_propertylist_append_double(aHeader, QC_ALIGN_DRA_MEAN, qcmean);
1355  cpl_propertylist_append_double(aHeader, QC_ALIGN_DRA_STDEV, qcsdev);
1356 
1357  qcmin = cpl_table_get_column_min(offsets, MUSE_OFFSETS_DDEC) * 3600.;
1358  qcmax = cpl_table_get_column_max(offsets, MUSE_OFFSETS_DDEC) * 3600.;
1359  qcmean = cpl_table_get_column_mean(offsets, MUSE_OFFSETS_DDEC) * 3600.;
1360  qcsdev = cpl_table_get_column_stdev(offsets, MUSE_OFFSETS_DDEC) * 3600.;
1361 
1362  cpl_propertylist_append_double(aHeader, QC_ALIGN_DDEC_MIN, qcmin);
1363  cpl_propertylist_append_double(aHeader, QC_ALIGN_DDEC_MAX, qcmax);
1364  cpl_propertylist_append_double(aHeader, QC_ALIGN_DDEC_MEAN, qcmean);
1365  cpl_propertylist_append_double(aHeader, QC_ALIGN_DDEC_STDEV, qcsdev);
1366 
1367  return offsets;
1368 }
1369 
1370 /*----------------------------------------------------------------------------*/
1379 /*----------------------------------------------------------------------------*/
1380 int
1381 muse_exp_align_compute(muse_processing *aProcessing,
1382  muse_exp_align_params_t *aParams)
1383 {
1384 
1385  cpl_size nexposures = cpl_frameset_count_tags(aProcessing->inframes,
1386  MUSE_TAG_IMAGE_FOV);
1387  if (nexposures < 2) {
1388  cpl_msg_error(__func__, "This recipe requires at least 2 input FOV "
1389  "images!");
1390  return -1;
1391  }
1392 
1393  /* Validate input parameters */
1394  cpl_msg_debug(__func__, "Validating source detection parameters");
1395  if (!muse_align_check_detection_params(aParams)) {
1396  cpl_error_set_where(__func__);
1397  return -1;
1398  }
1399 
1400  /* Load and validate all fov input images */
1401  cpl_msg_info(__func__, "Loading FOV images");
1402 
1403  muse_imagelist *fovimages = muse_processing_fov_load_all(aProcessing);
1404  if (!fovimages) {
1405  cpl_msg_error(__func__, "Could not load FOV input images!");
1406  return -1;
1407  }
1408  nexposures = (cpl_size)muse_imagelist_get_size(fovimages);
1409 
1410  cpl_size iexposure;
1411  for (iexposure = 0; iexposure < nexposures; ++iexposure) {
1412  const muse_image *fov = muse_imagelist_get(fovimages, iexposure);
1413  cpl_errorstate es = cpl_errorstate_get();
1414  double mjdobs = muse_pfits_get_mjdobs(fov->header);
1415 
1416  if (!cpl_errorstate_is_equal(es) &&
1417  cpl_error_get_code() == CPL_ERROR_DATA_NOT_FOUND) {
1418  cpl_msg_error(__func__, "Missing \"MJD-OBS\" in FOV image %"
1419  CPL_SIZE_FORMAT, iexposure + 1);
1420  muse_imagelist_delete(fovimages);
1421  return -1;
1422  } else {
1423  cpl_size jexposure;
1424  for (jexposure = iexposure + 1; jexposure < nexposures; ++jexposure) {
1425  const muse_image *_fov = muse_imagelist_get(fovimages, jexposure);
1426  double _mjdobs = muse_pfits_get_mjdobs(_fov->header);
1427  if (_mjdobs == mjdobs) {
1428  cpl_msg_warning(__func__, "Found FOV images with non-unique "
1429  "\"MJD-OBS\" value (image %" CPL_SIZE_FORMAT
1430  " and %" CPL_SIZE_FORMAT ")!", iexposure + 1,
1431  jexposure + 1);
1432  }
1433  }
1434  }
1435  } /* for iexposure */
1436 
1437  cpl_msg_info(__func__, "Computing pointing corrections for %" CPL_SIZE_FORMAT
1438  " FOV images", nexposures);
1439 
1440  /* Cleanup NaN values which are possibly present in the input images to *
1441  * flag bad pixel. Replace them with an appropriate numerical value */
1442  for (iexposure = 0; iexposure < nexposures; ++iexposure) {
1443  muse_image *image = muse_imagelist_get(fovimages, iexposure);
1444  muse_utils_replace_nan(image, 0.);
1445  }
1446 
1447  /* For each FOV input image run the source detection and collect *
1448  * the created source lists */
1449  int nobjects[2] = {aParams->srcmin, aParams->srcmax};
1450 
1451  double roundness[2] = {aParams->roundmin, aParams->roundmax};
1452  double sharpness[2] = {aParams->sharpmin, aParams->sharpmax};
1453 
1454  cpl_table **srclists = cpl_calloc(nexposures, sizeof *srclists);
1455 
1456  for (iexposure = 0; iexposure < nexposures; ++iexposure) {
1457  cpl_msg_debug(__func__, "Detecting sources on image %" CPL_SIZE_FORMAT
1458  " of %" CPL_SIZE_FORMAT, iexposure + 1, nexposures);
1459 
1460  muse_image *fov = muse_imagelist_get(fovimages, iexposure);
1461  double threshold = aParams->threshold;
1462  cpl_boolean iterate = CPL_TRUE;
1463  int iteration = 0;
1464  cpl_size ndetections = 0;
1465  cpl_table *detections;
1466  cpl_image *image = cpl_image_cast(fov->data, CPL_TYPE_DOUBLE);
1467 
1468  /* If the image does not have an even number of pixels along the *
1469  * any axis discard the last column and ro respectively. */
1470 
1471  /* XXX: The size limit along the x- and y-axis corresponds to *
1472  * minimum convolution box size of muse_find_stars(). This *
1473  * check should better be done in muse_find_stars() itself! */
1474  cpl_size nx = cpl_image_get_size_x(image);
1475  cpl_size ny = cpl_image_get_size_y(image);
1476 
1477  if ((nx < 6) || (ny < 6)) {
1478  cpl_image_delete(image);
1479  muse_vfree((void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
1480  muse_imagelist_delete(fovimages);
1481  cpl_msg_error(__func__, "Image %" CPL_SIZE_FORMAT " of %"
1482  CPL_SIZE_FORMAT "is too small!", iexposure + 1, nexposures);
1483  return -1;
1484  }
1485 
1486  if ((nx % 2 != 0) || (ny % 2 != 0)) {
1487  cpl_msg_debug(__func__, "Trimming image %" CPL_SIZE_FORMAT " to an "
1488  "even number of pixels along the axes.", iexposure + 1);
1489  if (nx % 2 != 0) {
1490  --nx;
1491  }
1492  if (ny % 2 != 0) {
1493  --ny;
1494  }
1495  cpl_image *_image = cpl_image_extract(image, 1, 1, nx, ny);
1496  if (!_image) {
1497  cpl_image_delete(image);
1498  muse_vfree((void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
1499  muse_imagelist_delete(fovimages);
1500  cpl_msg_error(__func__, "Trimming image %" CPL_SIZE_FORMAT " failed",
1501  iexposure + 1);
1502  return -1;
1503  }
1504  cpl_image_delete(image);
1505  image = _image;
1506  } /* if odd image dimension(s) */
1507 
1508  while (iterate) {
1509  cpl_errorstate status = cpl_errorstate_get();
1510 
1511  detections = muse_find_stars(image, threshold, aParams->fwhm,
1512  roundness, sharpness);
1513  if (!detections) {
1514  if (cpl_error_get_code() == CPL_ERROR_DATA_NOT_FOUND) {
1515  cpl_errorstate_set(status);
1516  ndetections = 0;
1517  } else {
1518  cpl_image_delete(image);
1519  muse_vfree((void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
1520  muse_imagelist_delete(fovimages);
1521  cpl_msg_error(__func__, "Source detection failed on image %"
1522  CPL_SIZE_FORMAT " of %" CPL_SIZE_FORMAT, iexposure + 1,
1523  nexposures);
1524  return -1;
1525  }
1526  } else {
1527  ndetections = cpl_table_get_nrow(detections);
1528  }
1529 
1530  if ((ndetections > 0) &&
1531  ((ndetections >= nobjects[0]) && (ndetections <= nobjects[1]))) {
1532  iterate = CPL_FALSE;
1533  } else {
1534  cpl_table_delete(detections);
1535  detections = NULL;
1536  if (++iteration < aParams->iterations) {
1537  if (ndetections > nobjects[1]) {
1538  threshold += aParams->step;
1539  } else {
1540  threshold -= aParams->step;
1541  if (threshold <= 0.) {
1542  iterate = CPL_FALSE;
1543  }
1544  }
1545  } else {
1546  iterate = CPL_FALSE;
1547  }
1548  }
1549  } /* while iterate */
1550  cpl_image_delete(image);
1551 
1552  if (!detections) {
1553  muse_vfree((void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
1554  muse_imagelist_delete(fovimages);
1555  if (ndetections == 0) {
1556  cpl_msg_error(__func__, "No point-sources were detected in image %"
1557  CPL_SIZE_FORMAT, iexposure + 1);
1558  } else if (ndetections < nobjects[0]) {
1559  cpl_msg_error(__func__, "Too few (%" CPL_SIZE_FORMAT ") point-sources "
1560  "were detected in image %" CPL_SIZE_FORMAT, ndetections,
1561  iexposure + 1);
1562  } else if (ndetections > nobjects[1]) {
1563  cpl_msg_error(__func__, "Too many (%" CPL_SIZE_FORMAT ") point-sources "
1564  "were detected in image %" CPL_SIZE_FORMAT, ndetections,
1565  iexposure + 1);
1566  } else {
1567  cpl_msg_error(__func__, "Invalid number (%" CPL_SIZE_FORMAT ") of "
1568  "point-sources were detected in image %" CPL_SIZE_FORMAT,
1569  ndetections, iexposure + 1);
1570  }
1571  return -1;
1572  } /* if no detections */
1573 
1574  cpl_msg_info(__func__, "%" CPL_SIZE_FORMAT " point-sources detected in "
1575  "image %" CPL_SIZE_FORMAT, ndetections, iexposure + 1);
1576 
1577  /* Compute right ascension and declination from the pixel coordinates of *
1578  * each detected point-source using the WCS information of each exposure */
1579  if (muse_align_celestial_from_pixel(detections, fov->header)) {
1580  muse_vfree((void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
1581  muse_imagelist_delete(fovimages);
1582  cpl_msg_error(__func__, "Computing WCS coordinates from pixel positions "
1583  "failed for image %" CPL_SIZE_FORMAT " of %"
1584  CPL_SIZE_FORMAT, iexposure + 1, nexposures);
1585  return -1;
1586  }
1587  srclists[iexposure] = detections;
1588  } /* for iexposure */
1589 
1590 
1591  /* Compute exposure offsets */
1592  cpl_array *radii = muse_align_parse_valuelist(aParams->rsearch);
1593  if (!radii) {
1594  muse_vfree((void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
1595  muse_imagelist_delete(fovimages);
1596  cpl_msg_error(__func__, "Cannot compute field offsets! No valid search "
1597  "radius was given!");
1598  return -1;
1599  }
1600 
1601  cpl_msg_info(__func__, "Calculating field offsets...");
1602  cpl_propertylist *header = cpl_propertylist_new();
1603  cpl_table *offsets = muse_align_compute_field_offsets(fovimages, srclists, radii,
1604  aParams->nbins, aParams->weight,
1605  header);
1606  cpl_array_delete(radii);
1607  muse_vfree((void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
1608  muse_imagelist_delete(fovimages);
1609  if (!offsets) {
1610  cpl_msg_error(__func__, "Could not compute FOV offsets for %"
1611  CPL_SIZE_FORMAT " images", nexposures);
1612  cpl_propertylist_delete(header);
1613  return -1;
1614  }
1615 
1616  /* Write offset table product to file */
1617  cpl_error_code rc = muse_processing_save_table(aProcessing, -1, offsets,
1618  header, MUSE_TAG_OFFSET_LIST,
1620  cpl_propertylist_delete(header);
1621  cpl_table_delete(offsets);
1622 
1623  if (rc != CPL_ERROR_NONE) {
1624  return -1;
1625  }
1626 
1627  return 0;
1628 } /* muse_exp_align_compute() */
int srcmin
Minimum number of sources which must be found.
Structure definition for a collection of muse_images.
int weight
Use weighting.
void muse_image_delete(muse_image *aImage)
Deallocate memory associated to a muse_image object.
Definition: muse_image.c:85
const char * muse_pfits_get_extname(const cpl_propertylist *aHeaders)
find out the extension name
Definition: muse_pfits.c:205
double fwhm
FWHM in pixels of the convolution filter.
A structure containing a spatial two-axis WCS.
Definition: muse_wcs.h:105
cpl_size muse_pfits_get_naxis(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the size of a given axis
Definition: muse_pfits.c:242
cpl_image * data
the data extension
Definition: muse_image.h:46
double roundmax
Upper limit of the allowed point-source roundness.
cpl_image * stat
the statistics extension
Definition: muse_image.h:64
const char * muse_pfits_get_dateobs(const cpl_propertylist *aHeaders)
find out the date of observations
Definition: muse_pfits.c:364
const muse_cpltable_def muse_offset_list_def[]
void muse_imagelist_delete(muse_imagelist *aList)
Free the memory of the MUSE image list.
muse_image * muse_fov_load(const char *aFilename)
Load a FOV image into a MUSE image.
Definition: muse_utils.c:2864
const char * muse_pfits_get_bunit(const cpl_propertylist *aHeaders)
find out the unit string
Definition: muse_pfits.c:223
Structure definition of MUSE three extension FITS file.
Definition: muse_image.h:40
const char * muse_pfits_get_cunit(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS axis unit string
Definition: muse_pfits.c:491
cpl_propertylist * header
the FITS header
Definition: muse_image.h:72
cpl_boolean muse_processing_check_intags(muse_processing *aProcessing, const char *aTag, int aNChars)
Check that a tag is part of the input tags of a processing structure.
double step
Increment/decrement of the threshold value in subsequent iterations.
unsigned int muse_imagelist_get_size(muse_imagelist *aList)
Return the number of stored images.
cpl_matrix * muse_cplmatrix_extract_selected(const cpl_matrix *aMatrix, const cpl_array *aIndices)
Extract the elements given by the index array.
cpl_image * dq
the data quality extension
Definition: muse_image.h:56
cpl_table * muse_cpltable_new(const muse_cpltable_def *aDef, cpl_size aLength)
Create an empty table according to the specified definition.
muse_wcs * muse_wcs_new(cpl_propertylist *aHeader)
Create a new WCS structure from a given FITS header.
Definition: muse_wcs.c:1870
cpl_array * muse_cplarray_new_from_delimited_string(const char *aString, const char *aDelim)
Convert a delimited string into an array of strings.
double sharpmax
Upper limit of the allowed point-source sharpness.
Structure to hold the parameters of the muse_exp_align recipe.
int srcmax
Maximum number of sources which may be found.
#define MUSE_WCS_KEYS
regular expression for WCS properties
Definition: muse_wcs.h:48
muse_image * muse_imagelist_get(muse_imagelist *aList, unsigned int aIdx)
Get the muse_image of given list index.
double threshold
Initial intensity threshold for detecting point sources in sigma above background RMS...
static void muse_wcs_celestial_from_pixel_fast(muse_wcs *aWCS, double aX, double aY, double *aRA, double *aDEC)
Convert from pixel coordinates to celestial spherical coordinates.
Definition: muse_wcs.h:138
int nbins
Number of bins to use for 2D histogram on the first iteration of the offset computation.
void muse_processing_append_used(muse_processing *aProcessing, cpl_frame *aFrame, cpl_frame_group aGroup, int aDuplicate)
Add a frame to the set of used frames.
const char * rsearch
Search radius (in arcsec) for each iteration of the offset computation.
double muse_pfits_get_mjdobs(const cpl_propertylist *aHeaders)
find out the Julian Date of the observation
Definition: muse_pfits.c:346
double roundmin
Lower limit of the allowed point-source roundness.
muse_imagelist * muse_imagelist_new(void)
Create a new (empty) MUSE image list.
cpl_matrix * muse_cplmatrix_multiply_create(const cpl_matrix *aMatrix1, const cpl_matrix *aMatrix2)
Compute the element-wise product of two matrices.
cpl_frameset * inframes
cpl_array * muse_cplmatrix_where(const cpl_matrix *aMatrix, double aValue, muse_cplmatrix_element_compare_func aCondition)
Select matrix elements according to a condition.
cpl_table * muse_find_stars(const cpl_image *aImage, double aHmin, double aFwhm, const double *aRoundLimits, const double *aSharpLimits)
Find positive brightness perturbations (i.e stars) in an image.
cpl_boolean iscelsph
Definition: muse_wcs.h:112
cpl_error_code muse_processing_save_table(muse_processing *aProcessing, int aIFU, void *aTable, cpl_propertylist *aHeader, const char *aTag, muse_table_type aType)
Save a computed table to disk.
muse_image * muse_image_new(void)
Allocate memory for a new muse_image object.
Definition: muse_image.c:66
cpl_error_code muse_imagelist_set(muse_imagelist *aList, muse_image *aImage, unsigned int aIdx)
Set the muse_image of given list index.
cpl_error_code muse_image_dq_to_nan(muse_image *aImage)
Convert pixels flagged in the DQ extension to NANs in DATA (and STAT, if present).
Definition: muse_image.c:904
double sharpmin
Lower limit of the allowed point-source sharpness.