UVES Pipeline Reference Manual  5.4.0
uves_wavecal_search.c
1 /* *
2  * This file is part of the ESO UVES Pipeline *
3  * Copyright (C) 2004,2005 European Southern Observatory *
4  * *
5  * This library is free software; you can redistribute it and/or modify *
6  * it under the terms of the GNU General Public License as published by *
7  * the Free Software Foundation; either version 2 of the License, or *
8  * (at your option) any later version. *
9  * *
10  * This program is distributed in the hope that it will be useful, *
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of *
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
13  * GNU General Public License for more details. *
14  * *
15  * You should have received a copy of the GNU General Public License *
16  * along with this program; if not, write to the Free Software *
17  * Foundation, 51 Franklin St, Fifth Floor, Boston, MA 02111-1307 USA *
18  * */
19 
20 /*
21  * $Author: amodigli $
22  * $Date: 2013-10-09 12:11:14 $
23  * $Revision: 1.34 $
24  * $Name: not supported by cvs2svn $
25  * $Log: not supported by cvs2svn $
26  * Revision 1.33 2012/03/02 16:43:12 amodigli
27  * fixed warning related to upgrade to CPL6
28  *
29  * Revision 1.32 2011/12/08 13:59:43 amodigli
30  * Fox warnings with CPL6
31  *
32  * Revision 1.31 2011/03/25 07:44:19 amodigli
33  * cleaned output
34  *
35  * Revision 1.30 2011/03/23 12:27:31 amodigli
36  * changed QC key as user likes
37  *
38  * Revision 1.29 2011/03/23 10:08:31 amodigli
39  * added QC to better characterize wave accuracy
40  *
41  * Revision 1.28 2010/09/24 09:32:10 amodigli
42  * put back QFITS dependency to fix problem spot by NRI on FIBER mode (with MIDAS calibs) data
43  *
44  * Revision 1.26 2007/08/21 13:08:26 jmlarsen
45  * Removed irplib_access module, largely deprecated by CPL-4
46  *
47  * Revision 1.25 2007/06/06 08:17:34 amodigli
48  * replace tab with 4 spaces
49  *
50  * Revision 1.24 2007/05/02 13:20:01 jmlarsen
51  * Take error bars into account in line searching if arclamp was flat-fielded
52  *
53  * Revision 1.23 2007/04/24 12:50:29 jmlarsen
54  * Replaced cpl_propertylist -> uves_propertylist which is much faster
55  *
56  * Revision 1.22 2007/04/20 14:46:45 jmlarsen
57  * Added commented out code
58  *
59  * Revision 1.21 2007/03/05 10:25:08 jmlarsen
60  * Include slope in Gaussian fit
61  *
62  * Revision 1.20 2007/02/23 13:33:38 jmlarsen
63  * Added code to test unweighted fitting
64  *
65  * Revision 1.19 2007/02/22 15:38:26 jmlarsen
66  * Use linear background term in Gaussian fit
67  *
68  * Revision 1.18 2006/11/15 15:02:15 jmlarsen
69  * Implemented const safe workarounds for CPL functions
70  *
71  * Revision 1.16 2006/11/15 14:04:08 jmlarsen
72  * Removed non-const version of parameterlist_get_first/last/next which is already
73  * in CPL, added const-safe wrapper, unwrapper and deallocator functions
74  *
75  * Revision 1.15 2006/11/06 15:19:42 jmlarsen
76  * Removed unused include directives
77  *
78  * Revision 1.14 2006/08/18 13:51:01 jmlarsen
79  * Moved one message from info to debug level
80  *
81  * Revision 1.13 2006/08/17 13:56:53 jmlarsen
82  * Reduced max line length
83  *
84  * Revision 1.12 2006/08/17 09:18:47 jmlarsen
85  * Removed CPL2 code
86  *
87  * Revision 1.11 2006/08/11 14:38:24 jmlarsen
88  * Minor text output change
89  *
90  * Revision 1.10 2006/08/11 09:01:17 jmlarsen
91  * Set unextracted bins to zero flux rather than marking them as bad
92  *
93  * Revision 1.9 2006/07/14 12:45:58 jmlarsen
94  * Removed a few messages
95  *
96  * Revision 1.8 2006/07/03 13:29:30 jmlarsen
97  * Adapted to new 1d-fitting function interface
98  *
99  * Revision 1.7 2006/06/13 12:05:11 jmlarsen
100  * Shortened max line length
101  *
102  * Revision 1.6 2006/05/12 15:06:30 jmlarsen
103  * Killed code for method = gravity
104  *
105  * Revision 1.5 2006/04/24 09:34:26 jmlarsen
106  * Adapted to new interface of gaussian fitting routine
107  *
108  * Revision 1.4 2006/03/03 13:54:11 jmlarsen
109  * Changed syntax of check macro
110  *
111  * Revision 1.3 2006/02/15 13:19:15 jmlarsen
112  * Reduced source code max. line length
113  *
114  * Revision 1.2 2006/02/08 07:52:16 jmlarsen
115  * Added function returning library version
116  *
117  * Revision 1.34 2006/01/12 15:41:14 jmlarsen
118  * Moved gauss. fitting to irplib
119  *
120  * Revision 1.33 2005/12/20 08:11:44 jmlarsen
121  * Added CVS entry
122  *
123  */
124 
125 /*----------------------------------------------------------------------------*/
129 /*----------------------------------------------------------------------------*/
132 #ifdef HAVE_CONFIG_H
133 # include <config.h>
134 #endif
135 
136 #include <uves_wavecal_search.h>
137 #include <uves_utils.h>
138 #include <uves_utils_wrappers.h>
139 #include <uves_utils_cpl.h>
140 #include <uves_pfits.h>
141 #include <uves_dump.h>
142 #include <uves_error.h>
143 #include <uves_msg.h>
144 #include <uves_qclog.h>
145 
146 #include <cpl.h>
147 #include <float.h>
148 
149 #define FIT_SLOPE 1
150 #define WEIGHTED_FIT 1 /* Define to zero to get unweighted fit of emmision line
151  (like MIDAS) */
152 
153 static double
154 xcenter(const cpl_image *image, const cpl_image *noise, int xlo, int xhi, int row,
155  centering_method CENTERING_METHOD, int bin_disp,
156  double *sigma, double *intensity, double *dx0, double *slope, double *background);
157 
158 static cpl_error_code
159 detect_lines(const cpl_image *spectrum, const cpl_image *noise,
160  const uves_propertylist *spectrum_header,
161  bool flat_fielded,
162  int RANGE, double THRESHOLD, centering_method CENTERING_METHOD,
163  int bin_disp,
164  const polynomial *order_locations, cpl_image *arcframe,
165  cpl_table *linetable,
166  int *ndetected, int *nrows);
167 
168 /*----------------------------------------------------------------------------*/
198 /*----------------------------------------------------------------------------*/
199 cpl_table *
200 uves_wavecal_search(const cpl_image *spectrum, const cpl_image *noise,
201  const uves_propertylist *spectrum_header,
202  bool flat_fielded,
203  const polynomial *order_locations, cpl_image *arcframe,
204  int RANGE, int MINLINES, int MAXLINES,
205  centering_method CENTERING_METHOD,int bin_disp,
206  const int trace, const int window, cpl_table* qclog)
207 {
208  cpl_table *linetable = NULL; /* Result */
209 
210  int nx, ny, norders; /* Dimensions of raw image, number of orders */
211  double threshold_low; /* Threshold limits used for binary search */
212  double threshold_high;
213  double threshold = 0; /* Current threshold */
214  int lines_in_table; /* Number of lines in line table */
215  int lines_detected; /* Number of lines actually found */
216  bool max_thresh_found = false; /* Is 'threshold_high' large enough? */
217 
218 
219  passure( spectrum != NULL, "Null input spectrum");
220  passure( order_locations != NULL, "Null polynomial");
221  passure( arcframe != NULL, "Null raw image");
222 
223  if (flat_fielded) {
224  assure( cpl_image_get_type(spectrum) == CPL_TYPE_DOUBLE,
225  CPL_ERROR_TYPE_MISMATCH,
226  "Spectrum image type is %s, must be double",
227  uves_tostring_cpl_type(cpl_image_get_type(spectrum)));
228  }
229 
230  check(( nx = cpl_image_get_size_x(spectrum),
231  norders = cpl_image_get_size_y(spectrum)), "Error reading input spectrum");
232  check( ny = cpl_image_get_size_y(arcframe), "Error reading input image");
233  assure(nx == cpl_image_get_size_x(arcframe), CPL_ERROR_INCOMPATIBLE_INPUT,
234  "Spectrum and image widths are different (%d and %" CPL_SIZE_FORMAT ")",
235  nx, cpl_image_get_size_x(arcframe));
236 
237  assure( MINLINES <= MAXLINES, CPL_ERROR_ILLEGAL_INPUT,
238  "minlines=%d maxlines=%d", MINLINES, MAXLINES );
239 
240  /* Initialize result line table */
241  check(( linetable = cpl_table_new(MAXLINES),
242  cpl_table_new_column(linetable, "X" , CPL_TYPE_DOUBLE),
243  cpl_table_new_column(linetable, "dX" , CPL_TYPE_DOUBLE),
244  cpl_table_new_column(linetable, "Xwidth", CPL_TYPE_DOUBLE),
245  cpl_table_new_column(linetable, "Y" , CPL_TYPE_INT),
246  cpl_table_new_column(linetable, "Peak" , CPL_TYPE_DOUBLE),
247  cpl_table_new_column(linetable, "Background" , CPL_TYPE_DOUBLE),
248  cpl_table_new_column(linetable, "Slope" , CPL_TYPE_DOUBLE)),
249  "Could not create line table");
250 
251  uves_msg("Searching for emission lines");
252 
253  threshold_low = 0.0;
254 
255  /* This start (guess) value is doubled until too few lines are detected */
256  if (flat_fielded) {
257  threshold_high = 10.0; /* dimensionless, number of stdevs above continuum */
258  }
259  else {
260  threshold_high = cpl_image_get_mean(spectrum);
261 
262  assure( threshold_high > 0, CPL_ERROR_ILLEGAL_INPUT,
263  "Spectrum median flux is %e. Must be positive",
264  cpl_image_get_median(spectrum));
265  }
266 
267  max_thresh_found = false;
268 
269  /* Detect lines and adjust threshold
270  until MINLINES <= lines_detected <= MAXLINES */
271  lines_detected = 0;
272 
273 
274  char qc_key[40];
275 
276  int kk=0;
277  while( (lines_detected < MINLINES || MAXLINES < lines_detected) &&
278  fabs(threshold_low - threshold_high) > DBL_EPSILON )
279  {
280  kk++;
281  threshold = (threshold_low + threshold_high)/2.0;
282 
283  check( detect_lines(spectrum, noise, spectrum_header,
284  flat_fielded,
285  RANGE, threshold, CENTERING_METHOD,
286  bin_disp,
287  order_locations,
288  NULL,
289  linetable,
290  &lines_detected,
291  &lines_in_table),
292  "Could not search for emission lines");
293 
294  /* Update threshold */
295  /* 'threshold_high' is doubled until the threshold is too high.
296  Then a binary search is performed. */
297  if (lines_detected < MINLINES)
298  {
299  max_thresh_found = true;
300  threshold_high = threshold;
301  }
302  else if (MAXLINES < lines_detected)
303  {
304  if (!max_thresh_found)
305  {
306  threshold_high *= 2;
307  }
308  else
309  {
310  threshold_low = threshold;
311  }
312  }
313  sprintf(qc_key,"QC TRACE%d WIN%d NLINDET%d",trace,window,kk);
314  uves_msg_debug("ThAr lamp on trace %d window %d detected lines %d",
315  trace,window,lines_detected);
316  ck0_nomsg(uves_qclog_add_int(qclog,qc_key,lines_detected,
317  "ThAr lamp detected lines","%d"));
318 
319 
320 
321 
322  } /* end while loop */
323 
324  sprintf(qc_key,"QC TRACE%d WIN%d NLINDET NITERS",trace,window);
325  ck0_nomsg(uves_qclog_add_int(qclog,qc_key,kk,"Number of iterations","%d") );
326  assure( MINLINES <= lines_detected && lines_detected <= MAXLINES,
327  CPL_ERROR_CONTINUE,
328  "Could not detect between %d and %d lines. Try to increase search range",
329  MINLINES, MAXLINES);
330 
331  /* Draw detections on input image */
332  check( detect_lines(spectrum, noise, spectrum_header,
333  flat_fielded,
334  RANGE, threshold, CENTERING_METHOD,
335  bin_disp,
336  order_locations,
337  arcframe,
338  linetable,
339  &lines_detected,
340  &lines_in_table),
341  "Could not search for emission lines");
342 
343  /* Remove the last part of the line table (garbage) */
344  check( cpl_table_set_size(linetable, lines_in_table),
345  "Could not resize line table");
346 
347  uves_sort_table_1(linetable, "X", false);
348 
349  cleanup:
350 #if 0 /* if flat-fielded */
351  uves_free_image(&temp);
352 #endif
353  if (cpl_error_get_code() != CPL_ERROR_NONE)
354  {
355  uves_free_table(&linetable);
356  }
357  else
358  {
359  /* Returned is... */
360  passure( cpl_table_get_ncol(linetable) == 7, "%" CPL_SIZE_FORMAT "",
361  cpl_table_get_ncol(linetable));
362  passure( cpl_table_has_column(linetable, "X" ), " ");
363  passure( cpl_table_has_column(linetable, "dX" ), " ");
364  passure( cpl_table_has_column(linetable, "Xwidth"), " ");
365  passure( cpl_table_has_column(linetable, "Y" ), " ");
366  passure( cpl_table_has_column(linetable, "Peak" ), " ");
367  passure( cpl_table_has_column(linetable, "Background" ), " ");
368  passure( cpl_table_has_column(linetable, "Slope" ), " ");
369 
370  }
371  return linetable;
372 }
373 
374 /*----------------------------------------------------------------------------*/
425 /*----------------------------------------------------------------------------*/
426 static cpl_error_code
427 detect_lines(const cpl_image *spectrum, const cpl_image *noise,
428  const uves_propertylist *spectrum_header,
429  bool flat_fielded,
430  int RANGE, double THRESHOLD, centering_method CENTERING_METHOD,
431  int bin_disp,
432  const polynomial *order_locations, cpl_image *arcframe,
433  cpl_table *linetable,
434  int *ndetected, int *nrows)
435 {
436  int norders; /* Number of orders */
437  int minorder; /* Relative order number of first row in spectrum image */
438  int MAXLINES; /* Number of rows in line table (max no. of
439  lines to search for) */
440  int nx; /* Width of spectrum (and raw image) */
441  int x, order; /* 'order' always counts from 1 */
442 
443  const double *spectrum_data;
444  const double *noise_data;
445 
446  /* Check input */
447  passure( spectrum != NULL, " ");
448  passure( noise != NULL, " ");
449  passure( spectrum_header != NULL, " ");
450  nx = cpl_image_get_size_x(spectrum);
451  norders = cpl_image_get_size_y(spectrum);
452 
453  /* For efficiency reasons, get direct pointer to buffer,
454  support only CPL_TYPE_DOUBLE */
455  assure( cpl_image_get_type(spectrum) == CPL_TYPE_DOUBLE,
456  CPL_ERROR_UNSUPPORTED_MODE,
457  "Image type must be double. It is %s",
458  uves_tostring_cpl_type(cpl_image_get_type(spectrum)));
459 
460  spectrum_data = cpl_image_get_data_double_const(spectrum);
461  noise_data = cpl_image_get_data_double_const(noise);
462 
463  passure( RANGE > 0, "%d", RANGE);
464 
465  if (arcframe != NULL)
466  {
467  passure( order_locations != NULL, " ");
468  passure( nx == cpl_image_get_size_x(arcframe),
469  "%d %" CPL_SIZE_FORMAT "", nx, cpl_image_get_size_x(arcframe));
470  }
471 
472  passure( linetable != NULL, " ");
473  MAXLINES = cpl_table_get_nrow(linetable);
474  passure( cpl_table_get_ncol(linetable) == 7, "%" CPL_SIZE_FORMAT "",
475  cpl_table_get_ncol(linetable));
476  passure( cpl_table_has_column(linetable, "X" ), " ");
477  passure( cpl_table_has_column(linetable, "dX" ), " ");
478  passure( cpl_table_has_column(linetable, "Xwidth"), " ");
479  passure( cpl_table_has_column(linetable, "Y" ), " ");
480  passure( cpl_table_has_column(linetable, "Peak" ), " ");
481  passure( cpl_table_has_column(linetable, "Background" ), " ");
482  passure( cpl_table_has_column(linetable, "Slope" ), " ");
483 
484  assure( THRESHOLD > 0, CPL_ERROR_ILLEGAL_INPUT, "Illegal threshold: %e",
485  THRESHOLD);
486 
487  check( minorder = uves_pfits_get_crval2(spectrum_header),
488  "Error reading order number of first row");
489 
490  *ndetected = 0; /* Counts the number of lines detected so far. */
491  *nrows = 0; /* A pointer to the first unused row in the
492  line table */
493 
494  for (order = minorder; order < minorder + norders; order++) {
495  int spectrum_row = order - minorder + 1;
496  int ndetected_order = 0;
497  for (x = 1; x <= nx; x++) {
498  double flux, dflux;
499  int peak_width = 0;
500  int xlo, xhi;
501  double local_median;
502 
503  /* Check if there is a peak and determine its position and width */
504  // flux = cpl_image_get(spectrum, x, spectrum_row, &pis_rejected);
505  flux = spectrum_data[(x-1) + (spectrum_row - 1) * nx];
506  dflux = noise_data [(x-1) + (spectrum_row - 1) * nx];
507 
508  xlo = uves_max_int(x - RANGE, 1);
509  xhi = uves_min_int(x + RANGE, nx);
510 
511  local_median = cpl_image_get_median_window(
512  spectrum,
513  uves_max_int(xlo, 1 ), spectrum_row,
514  uves_min_int(xhi, nx), spectrum_row);
515 
516  while(x <= nx &&
517  (
518  (!flat_fielded && flux - local_median > THRESHOLD)
519  ||
520  (flat_fielded && (flux - local_median) > THRESHOLD * dflux)
521  )
522  ) {
523 #if WANT_BIG_LOGFILE
524  uves_msg_debug("threshold = %f\tx = %d\tflux = %f\tmedian = %f",
525  THRESHOLD, x, flux, local_median);
526 #endif
527 
528  x += 1;
529  peak_width += 1;
530 
531  if (x <= nx) {
532  /* flux =
533  cpl_image_get(spectrum, x,
534  spectrum_row, &pis_rejected);
535  */
536  flux = spectrum_data[(x-1) + (spectrum_row - 1) * nx];
537  xlo = uves_max_int(x - RANGE, 1);
538  xhi = uves_min_int(x + RANGE, nx);
539  local_median = cpl_image_get_median_window(
540  spectrum,
541  uves_max_int(xlo, 1 ), spectrum_row,
542  uves_min_int(xhi, nx), spectrum_row);
543  }
544  }
545  /* x is now first position that is below (median + threshold) */
546 
547  if (peak_width > 0) {
548  double x_peak, dx = 0, sigma, slope, back;
549  check( x_peak = xcenter(spectrum, noise,
550  uves_max_int(1, x - peak_width),
551  /* First position above threshold */
552  uves_max_int(1, x - 1),
553  /* Last position above threshold */
554  spectrum_row,
555  CENTERING_METHOD,
556  bin_disp,
557  &sigma,
558  &flux,
559  &dx,
560  &slope,
561  &back),
562  "Could not locate peak center");
563 
564 #if WANT_BIG_LOGFILE
565  uves_msg_debug("(Order, x, flux) = (%d, %f, %f)",
566  order, x_peak, flux);
567 #endif
568  /* Add line to line table, but only if less
569  lines that MAXLINES have been detected */
570  if (*nrows < MAXLINES) {
571  check(( cpl_table_set_int (linetable, "Y" , *nrows, order),
572  cpl_table_set_double(linetable, "X" , *nrows, x_peak),
573  cpl_table_set_double(linetable, "dX" , *nrows, dx),
574  cpl_table_set_double(linetable, "Xwidth", *nrows, sigma),
575  cpl_table_set_double(linetable, "Peak" , *nrows, flux),
576  cpl_table_set_double(linetable, "Background" , *nrows, back),
577  cpl_table_set_double(linetable, "Slope" , *nrows, slope)),
578  "Could not update line table row %d", *nrows);
579  (*nrows)++;
580  }
581 
582  ndetected_order++;
583  (*ndetected)++;
584 
585  if (arcframe != NULL) {
586  int x1;
587  int pen = 0; /* Value to write */
588  int ny = cpl_image_get_size_y(arcframe);
589  /* We already know 'nx' from the width of the spectrum image */
590 
591  for (x1 = uves_max_int(
592  1 , uves_round_double(
593  x_peak - peak_width - 0*RANGE/2.0));
594  x1 <= uves_min_int(
595  nx, uves_round_double(
596  x_peak + peak_width + 0*RANGE/2.0));
597  x1++) {
598  check( cpl_image_set(
599  arcframe,
600  x1,
601  uves_min_int(
602  ny,
603  uves_max_int(
604  1,
606  order_locations, x1, order)
607  )),
608  pen),
609  "Error writing input image");
610  check( cpl_image_set(
611  arcframe,
612  uves_min_int(
613  nx,
614  uves_max_int((int) x_peak, 1)),
615  uves_min_int(
616  ny,
617  uves_max_int(
618  1,
620  order_locations, x1, order)
621  - 10)),
622  pen),
623  "Error writing input image");
624  }
625  }
626  } /* line found */
627  }/* for x */
628  if (arcframe != NULL) uves_msg_debug("Order #%d: %d lines detected",
629  order, ndetected_order);
630  }/* for order */
631 
632  /* Remove doublets */
633  {
634  int i;
635  int doublets_removed = 0;
636  for (i = 0; i+1 < *nrows; i++) {
637  if (fabs(cpl_table_get_double(linetable, "X", i , NULL) -
638  cpl_table_get_double(linetable, "X", i+1, NULL)) < 2.0)
639  {
640  /* If a doublet was found, delete it.
641  Make sure the table stays the same size
642  by adding two rows at the end. */
643 
644  check( cpl_table_erase_window(linetable, i, 2),
645  "Error removing rows");
646  *nrows -= 2;
647  *ndetected -= 2;
648 
649  check( cpl_table_set_size(linetable,
650  cpl_table_get_nrow(linetable) + 2),
651  "Could not resize line table");
652 
653  doublets_removed++;
654  }
655  }
656  if (doublets_removed > 0)
657  {
658  uves_msg_debug("%d doublet%s removed",
659  doublets_removed, doublets_removed > 1 ? "s" : "");
660  }
661  }
662 
663  uves_msg("Range = %d pixels; threshold = %.2f %s; %d lines detected",
664  RANGE, THRESHOLD, flat_fielded ? "stdev" : "ADU", *ndetected);
665 
666  cleanup:
667  return cpl_error_get_code();
668 }
669 
670 /*----------------------------------------------------------------------------*/
695 /*----------------------------------------------------------------------------*/
696 static double
697 xcenter(const cpl_image *image, const cpl_image *noise, int xlo, int xhi, int row,
698  centering_method CENTERING_METHOD, int bin_disp,
699  double *sigma, double *intensity, double *dx0, double *slope, double *background)
700 {
701  double x0; /* Result */
702  cpl_matrix *covariance = NULL;
703  const double *image_data;
704  bool converged;
705  int lo_r, hi_r;
706 
707  int nx = cpl_image_get_size_x(image);
708 
709  passure(cpl_image_get_type(image) == CPL_TYPE_DOUBLE, " ");
710 
711  image_data = cpl_image_get_data_double_const(image);
712 
713  /* Make sure fit window is 13-19 pixels
714  (7-9 pixels for binned CCD) */
715  lo_r = 6;
716  hi_r = 8;
717  if (bin_disp >= 2)
718  {
719  lo_r = 4;
720  hi_r = 5;
721  }
722 
723  {
724  int xm = (xlo+xhi)/2;
725 
726  xlo = uves_max_int(1, xm - lo_r);
727  xhi = uves_min_int(nx, xm + lo_r);
728  }
729 
730  /* Increase fit window (up to 19 pixels) */
731  do {
732  converged = true;
733  if (1 < xlo && 0 <
734  //cpl_image_get(image, xlo - 1, row, &pis_rejected) &&
735  //cpl_image_get(image, xlo - 1, row, &pis_rejected) <
736  //cpl_image_get(image, xlo , row, &pis_rejected) )
737  image_data[(xlo-1-1) + (row - 1) * nx] &&
738  image_data[(xlo-1-1) + (row - 1) * nx] <
739  image_data[(xlo -1) + (row - 1) * nx] )
740  {
741  converged = false;
742  xlo -= 1;
743  }
744 
745  if (xhi < nx && 0 <
746  //cpl_image_get(image, xhi + 1, row, &pis_rejected) &&
747  //cpl_image_get(image, xhi + 1, row, &pis_rejected) <
748  //cpl_image_get(image, xhi , row, &pis_rejected) )
749  image_data[(xhi+1-1) + (row - 1) * nx] &&
750  image_data[(xhi+1-1) + (row - 1) * nx] <
751  image_data[(xhi -1) + (row - 1) * nx] )
752  {
753  converged = false;
754  xhi += 1;
755  }
756 
757  if ((xhi-xlo+1) >= hi_r)
758  {
759  converged = true;
760  }
761 
762  } while (!converged);
763 
764  /* Get precise location */
765  if (CENTERING_METHOD == CENTERING_GAUSSIAN)
766  {
767 #if WEIGHTED_FIT
768  uves_fit_1d_image(image, noise, NULL,
769 #else /* Unweighted fit like MIDAS which gives larger sigma */
770  uves_fit_1d_image(image, NULL, NULL,
771 #endif
772  true, false, false,
773  xlo, xhi, row,
774  &x0, sigma, intensity, background, slope,
775 #if WEIGHTED_FIT
776  NULL, NULL, &covariance,
777 #else
778  NULL, NULL, NULL,
779 #endif
780 
781 #if FIT_SLOPE
783 #else
785  *slope = 0;
786 #endif
787 
788  /* The fitting routine sometimes (i.e. regularly) fails
789  * because of low statistics.
790  * Recover from specific fitting errors */
791  if (cpl_error_get_code() == CPL_ERROR_NONE)
792  {
793  /* Variance is guaranteed to be positive */
794 #if WEIGHTED_FIT
795  *dx0 = sqrt(cpl_matrix_get(covariance, 0, 0));
796 #else
797  *dx0 = *sigma / sqrt(*intensity);
798 #endif
799 
800 #if WANT_BIG_LOGFILE
801  uves_msg_debug("Gaussian fit succeeded at (x, row, N) = (%f, %d, %d)",
802  x0, row, xhi-xlo+1);
803 #endif
804  }
805  else if (cpl_error_get_code() == CPL_ERROR_CONTINUE)
806  {
807  /* Fitting failed */
809 #if WANT_BIG_LOGFILE
810  uves_msg_debug("Gaussian fit failed at (x, row, N) ="
811  " (%f, %d, %d), using centroid",
812  x0, row, xhi-xlo+1);
813 #endif
814  *dx0 = *sigma / sqrt(*intensity);
815  }
816  else if (cpl_error_get_code() == CPL_ERROR_SINGULAR_MATRIX)
817  {
819 
820  /* Fitting succeeded but covariance computation failed */
821  uves_msg_debug("Covariance matrix computation failed");
822  *dx0 = *sigma / sqrt(*intensity);
823  }
824 
825  assure(cpl_error_get_code() == CPL_ERROR_NONE, cpl_error_get_code(),
826  "Gaussian fitting failed");
827 
828 #if WANT_BIG_LOGFILE
829  uves_msg_debug("Fit = (x0=%f, sigma=%f, norm=%f, backg=%f, N=%d)",
830  x0,
831  *sigma,
832  *intensity,
833  background,
834  xhi - xlo + 1);
835 #endif
836 
837  /* 'intensity' is the norm (area) of the gaussian fit.
838  But we need to return the height above zero
839  (for MIDAS compatibility).
840  height = f(x0) = background + norm/sqrt(2pi sigma^2)
841  */
842 
843  *intensity = *background + (*intensity)/(sqrt(2*M_PI) * (*sigma));
844 
845  }
846  else /* if (CENTERING_METHOD == CENTERING_GRAVITY) */
847  {
848  assure (false, CPL_ERROR_UNSUPPORTED_MODE,
849  "Centering method (no. %d) is unsupported",
850  CENTERING_METHOD);
851  }
852 
853  cleanup:
854  uves_free_matrix(&covariance);
855  return x0;
856 }
857