UVES Pipeline Reference Manual  5.4.6
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  cpl_table_set_column_unit(linetable,"X","pix");
251  cpl_table_set_column_unit(linetable,"dX","pix");
252  cpl_table_set_column_unit(linetable,"Xwidth","pix");
253  cpl_table_set_column_unit(linetable,"Y","pix");
254  cpl_table_set_column_unit(linetable,"Peak","ADU");
255  cpl_table_set_column_unit(linetable,"Background","ADU");
256  cpl_table_set_column_unit(linetable,"Slope"," ");
257  uves_msg("Searching for emission lines");
258 
259  threshold_low = 0.0;
260 
261  /* This start (guess) value is doubled until too few lines are detected */
262  if (flat_fielded) {
263  threshold_high = 10.0; /* dimensionless, number of stdevs above continuum */
264  }
265  else {
266  threshold_high = cpl_image_get_mean(spectrum);
267 
268  assure( threshold_high > 0, CPL_ERROR_ILLEGAL_INPUT,
269  "Spectrum median flux is %e. Must be positive",
270  cpl_image_get_median(spectrum));
271  }
272 
273  max_thresh_found = false;
274 
275  /* Detect lines and adjust threshold
276  until MINLINES <= lines_detected <= MAXLINES */
277  lines_detected = 0;
278 
279 
280  char qc_key[40];
281 
282  int kk=0;
283  while( (lines_detected < MINLINES || MAXLINES < lines_detected) &&
284  fabs(threshold_low - threshold_high) > DBL_EPSILON )
285  {
286  kk++;
287  threshold = (threshold_low + threshold_high)/2.0;
288 
289  check( detect_lines(spectrum, noise, spectrum_header,
290  flat_fielded,
291  RANGE, threshold, CENTERING_METHOD,
292  bin_disp,
293  order_locations,
294  NULL,
295  linetable,
296  &lines_detected,
297  &lines_in_table),
298  "Could not search for emission lines");
299 
300  /* Update threshold */
301  /* 'threshold_high' is doubled until the threshold is too high.
302  Then a binary search is performed. */
303  if (lines_detected < MINLINES)
304  {
305  max_thresh_found = true;
306  threshold_high = threshold;
307  }
308  else if (MAXLINES < lines_detected)
309  {
310  if (!max_thresh_found)
311  {
312  threshold_high *= 2;
313  }
314  else
315  {
316  threshold_low = threshold;
317  }
318  }
319  sprintf(qc_key,"QC TRACE%d WIN%d NLINDET%d",trace,window,kk);
320  uves_msg_debug("ThAr lamp on trace %d window %d detected lines %d",
321  trace,window,lines_detected);
322  ck0_nomsg(uves_qclog_add_int(qclog,qc_key,lines_detected,
323  "ThAr lamp detected lines","%d"));
324 
325 
326 
327 
328  } /* end while loop */
329 
330  sprintf(qc_key,"QC TRACE%d WIN%d NLINDET NITERS",trace,window);
331  ck0_nomsg(uves_qclog_add_int(qclog,qc_key,kk,"Number of iterations","%d") );
332  assure( MINLINES <= lines_detected && lines_detected <= MAXLINES,
333  CPL_ERROR_CONTINUE,
334  "Could not detect between %d and %d lines. Try to increase search range",
335  MINLINES, MAXLINES);
336 
337  /* Draw detections on input image */
338  check( detect_lines(spectrum, noise, spectrum_header,
339  flat_fielded,
340  RANGE, threshold, CENTERING_METHOD,
341  bin_disp,
342  order_locations,
343  arcframe,
344  linetable,
345  &lines_detected,
346  &lines_in_table),
347  "Could not search for emission lines");
348 
349  /* Remove the last part of the line table (garbage) */
350  check( cpl_table_set_size(linetable, lines_in_table),
351  "Could not resize line table");
352 
353  uves_sort_table_1(linetable, "X", false);
354 
355  cleanup:
356 #if 0 /* if flat-fielded */
357  uves_free_image(&temp);
358 #endif
359  if (cpl_error_get_code() != CPL_ERROR_NONE)
360  {
361  uves_free_table(&linetable);
362  }
363  else
364  {
365  /* Returned is... */
366  passure( cpl_table_get_ncol(linetable) == 7, "%" CPL_SIZE_FORMAT "",
367  cpl_table_get_ncol(linetable));
368  passure( cpl_table_has_column(linetable, "X" ), " ");
369  passure( cpl_table_has_column(linetable, "dX" ), " ");
370  passure( cpl_table_has_column(linetable, "Xwidth"), " ");
371  passure( cpl_table_has_column(linetable, "Y" ), " ");
372  passure( cpl_table_has_column(linetable, "Peak" ), " ");
373  passure( cpl_table_has_column(linetable, "Background" ), " ");
374  passure( cpl_table_has_column(linetable, "Slope" ), " ");
375 
376  }
377  return linetable;
378 }
379 
380 /*----------------------------------------------------------------------------*/
431 /*----------------------------------------------------------------------------*/
432 static cpl_error_code
433 detect_lines(const cpl_image *spectrum, const cpl_image *noise,
434  const uves_propertylist *spectrum_header,
435  bool flat_fielded,
436  int RANGE, double THRESHOLD, centering_method CENTERING_METHOD,
437  int bin_disp,
438  const polynomial *order_locations, cpl_image *arcframe,
439  cpl_table *linetable,
440  int *ndetected, int *nrows)
441 {
442  int norders; /* Number of orders */
443  int minorder; /* Relative order number of first row in spectrum image */
444  int MAXLINES; /* Number of rows in line table (max no. of
445  lines to search for) */
446  int nx; /* Width of spectrum (and raw image) */
447  int x, order; /* 'order' always counts from 1 */
448 
449  const double *spectrum_data;
450  const double *noise_data;
451 
452  /* Check input */
453  passure( spectrum != NULL, " ");
454  passure( noise != NULL, " ");
455  passure( spectrum_header != NULL, " ");
456  nx = cpl_image_get_size_x(spectrum);
457  norders = cpl_image_get_size_y(spectrum);
458 
459  /* For efficiency reasons, get direct pointer to buffer,
460  support only CPL_TYPE_DOUBLE */
461  assure( cpl_image_get_type(spectrum) == CPL_TYPE_DOUBLE,
462  CPL_ERROR_UNSUPPORTED_MODE,
463  "Image type must be double. It is %s",
464  uves_tostring_cpl_type(cpl_image_get_type(spectrum)));
465 
466  spectrum_data = cpl_image_get_data_double_const(spectrum);
467  noise_data = cpl_image_get_data_double_const(noise);
468 
469  passure( RANGE > 0, "%d", RANGE);
470 
471  if (arcframe != NULL)
472  {
473  passure( order_locations != NULL, " ");
474  passure( nx == cpl_image_get_size_x(arcframe),
475  "%d %" CPL_SIZE_FORMAT "", nx, cpl_image_get_size_x(arcframe));
476  }
477 
478  passure( linetable != NULL, " ");
479  MAXLINES = cpl_table_get_nrow(linetable);
480  passure( cpl_table_get_ncol(linetable) == 7, "%" CPL_SIZE_FORMAT "",
481  cpl_table_get_ncol(linetable));
482  passure( cpl_table_has_column(linetable, "X" ), " ");
483  passure( cpl_table_has_column(linetable, "dX" ), " ");
484  passure( cpl_table_has_column(linetable, "Xwidth"), " ");
485  passure( cpl_table_has_column(linetable, "Y" ), " ");
486  passure( cpl_table_has_column(linetable, "Peak" ), " ");
487  passure( cpl_table_has_column(linetable, "Background" ), " ");
488  passure( cpl_table_has_column(linetable, "Slope" ), " ");
489 
490  assure( THRESHOLD > 0, CPL_ERROR_ILLEGAL_INPUT, "Illegal threshold: %e",
491  THRESHOLD);
492 
493  check( minorder = uves_pfits_get_crval2(spectrum_header),
494  "Error reading order number of first row");
495 
496  *ndetected = 0; /* Counts the number of lines detected so far. */
497  *nrows = 0; /* A pointer to the first unused row in the
498  line table */
499 
500  for (order = minorder; order < minorder + norders; order++) {
501  int spectrum_row = order - minorder + 1;
502  int ndetected_order = 0;
503  for (x = 1; x <= nx; x++) {
504  double flux, dflux;
505  int peak_width = 0;
506  int xlo, xhi;
507  double local_median;
508 
509  /* Check if there is a peak and determine its position and width */
510  // flux = cpl_image_get(spectrum, x, spectrum_row, &pis_rejected);
511  flux = spectrum_data[(x-1) + (spectrum_row - 1) * nx];
512  dflux = noise_data [(x-1) + (spectrum_row - 1) * nx];
513 
514  xlo = uves_max_int(x - RANGE, 1);
515  xhi = uves_min_int(x + RANGE, nx);
516 
517  local_median = cpl_image_get_median_window(
518  spectrum,
519  uves_max_int(xlo, 1 ), spectrum_row,
520  uves_min_int(xhi, nx), spectrum_row);
521 
522  while(x <= nx &&
523  (
524  (!flat_fielded && flux - local_median > THRESHOLD)
525  ||
526  (flat_fielded && (flux - local_median) > THRESHOLD * dflux)
527  )
528  ) {
529 #if WANT_BIG_LOGFILE
530  uves_msg_debug("threshold = %f\tx = %d\tflux = %f\tmedian = %f",
531  THRESHOLD, x, flux, local_median);
532 #endif
533 
534  x += 1;
535  peak_width += 1;
536 
537  if (x <= nx) {
538  /* flux =
539  cpl_image_get(spectrum, x,
540  spectrum_row, &pis_rejected);
541  */
542  flux = spectrum_data[(x-1) + (spectrum_row - 1) * nx];
543  xlo = uves_max_int(x - RANGE, 1);
544  xhi = uves_min_int(x + RANGE, nx);
545  local_median = cpl_image_get_median_window(
546  spectrum,
547  uves_max_int(xlo, 1 ), spectrum_row,
548  uves_min_int(xhi, nx), spectrum_row);
549  }
550  }
551  /* x is now first position that is below (median + threshold) */
552 
553  if (peak_width > 0) {
554  double x_peak, dx = 0, sigma, slope, back;
555  check( x_peak = xcenter(spectrum, noise,
556  uves_max_int(1, x - peak_width),
557  /* First position above threshold */
558  uves_max_int(1, x - 1),
559  /* Last position above threshold */
560  spectrum_row,
561  CENTERING_METHOD,
562  bin_disp,
563  &sigma,
564  &flux,
565  &dx,
566  &slope,
567  &back),
568  "Could not locate peak center");
569 
570 #if WANT_BIG_LOGFILE
571  uves_msg_debug("(Order, x, flux) = (%d, %f, %f)",
572  order, x_peak, flux);
573 #endif
574  /* Add line to line table, but only if less
575  lines that MAXLINES have been detected */
576  if (*nrows < MAXLINES) {
577  check(( cpl_table_set_int (linetable, "Y" , *nrows, order),
578  cpl_table_set_double(linetable, "X" , *nrows, x_peak),
579  cpl_table_set_double(linetable, "dX" , *nrows, dx),
580  cpl_table_set_double(linetable, "Xwidth", *nrows, sigma),
581  cpl_table_set_double(linetable, "Peak" , *nrows, flux),
582  cpl_table_set_double(linetable, "Background" , *nrows, back),
583  cpl_table_set_double(linetable, "Slope" , *nrows, slope)),
584  "Could not update line table row %d", *nrows);
585  (*nrows)++;
586  }
587 
588  ndetected_order++;
589  (*ndetected)++;
590 
591  if (arcframe != NULL) {
592  int x1;
593  int pen = 0; /* Value to write */
594  int ny = cpl_image_get_size_y(arcframe);
595  /* We already know 'nx' from the width of the spectrum image */
596 
597  for (x1 = uves_max_int(
598  1 , uves_round_double(
599  x_peak - peak_width - 0*RANGE/2.0));
600  x1 <= uves_min_int(
601  nx, uves_round_double(
602  x_peak + peak_width + 0*RANGE/2.0));
603  x1++) {
604  check( cpl_image_set(
605  arcframe,
606  x1,
607  uves_min_int(
608  ny,
609  uves_max_int(
610  1,
612  order_locations, x1, order)
613  )),
614  pen),
615  "Error writing input image");
616  check( cpl_image_set(
617  arcframe,
618  uves_min_int(
619  nx,
620  uves_max_int((int) x_peak, 1)),
621  uves_min_int(
622  ny,
623  uves_max_int(
624  1,
626  order_locations, x1, order)
627  - 10)),
628  pen),
629  "Error writing input image");
630  }
631  }
632  } /* line found */
633  }/* for x */
634  if (arcframe != NULL) uves_msg_debug("Order #%d: %d lines detected",
635  order, ndetected_order);
636  }/* for order */
637 
638  /* Remove doublets */
639  {
640  int i;
641  int doublets_removed = 0;
642  for (i = 0; i+1 < *nrows; i++) {
643  if (fabs(cpl_table_get_double(linetable, "X", i , NULL) -
644  cpl_table_get_double(linetable, "X", i+1, NULL)) < 2.0)
645  {
646  /* If a doublet was found, delete it.
647  Make sure the table stays the same size
648  by adding two rows at the end. */
649 
650  check( cpl_table_erase_window(linetable, i, 2),
651  "Error removing rows");
652  *nrows -= 2;
653  *ndetected -= 2;
654 
655  check( cpl_table_set_size(linetable,
656  cpl_table_get_nrow(linetable) + 2),
657  "Could not resize line table");
658 
659  doublets_removed++;
660  }
661  }
662  if (doublets_removed > 0)
663  {
664  uves_msg_debug("%d doublet%s removed",
665  doublets_removed, doublets_removed > 1 ? "s" : "");
666  }
667  }
668 
669  uves_msg("Range = %d pixels; threshold = %.2f %s; %d lines detected",
670  RANGE, THRESHOLD, flat_fielded ? "stdev" : "ADU", *ndetected);
671 
672  cleanup:
673  return cpl_error_get_code();
674 }
675 
676 /*----------------------------------------------------------------------------*/
701 /*----------------------------------------------------------------------------*/
702 static double
703 xcenter(const cpl_image *image, const cpl_image *noise, int xlo, int xhi, int row,
704  centering_method CENTERING_METHOD, int bin_disp,
705  double *sigma, double *intensity, double *dx0, double *slope, double *background)
706 {
707  double x0; /* Result */
708  cpl_matrix *covariance = NULL;
709  const double *image_data;
710  bool converged;
711  int lo_r, hi_r;
712 
713  int nx = cpl_image_get_size_x(image);
714 
715  passure(cpl_image_get_type(image) == CPL_TYPE_DOUBLE, " ");
716 
717  image_data = cpl_image_get_data_double_const(image);
718 
719  /* Make sure fit window is 13-19 pixels
720  (7-9 pixels for binned CCD) */
721  lo_r = 6;
722  hi_r = 8;
723  if (bin_disp >= 2)
724  {
725  lo_r = 4;
726  hi_r = 5;
727  }
728 
729  {
730  int xm = (xlo+xhi)/2;
731 
732  xlo = uves_max_int(1, xm - lo_r);
733  xhi = uves_min_int(nx, xm + lo_r);
734  }
735 
736  /* Increase fit window (up to 19 pixels) */
737  do {
738  converged = true;
739  if (1 < xlo && 0 <
740  //cpl_image_get(image, xlo - 1, row, &pis_rejected) &&
741  //cpl_image_get(image, xlo - 1, row, &pis_rejected) <
742  //cpl_image_get(image, xlo , row, &pis_rejected) )
743  image_data[(xlo-1-1) + (row - 1) * nx] &&
744  image_data[(xlo-1-1) + (row - 1) * nx] <
745  image_data[(xlo -1) + (row - 1) * nx] )
746  {
747  converged = false;
748  xlo -= 1;
749  }
750 
751  if (xhi < nx && 0 <
752  //cpl_image_get(image, xhi + 1, row, &pis_rejected) &&
753  //cpl_image_get(image, xhi + 1, row, &pis_rejected) <
754  //cpl_image_get(image, xhi , row, &pis_rejected) )
755  image_data[(xhi+1-1) + (row - 1) * nx] &&
756  image_data[(xhi+1-1) + (row - 1) * nx] <
757  image_data[(xhi -1) + (row - 1) * nx] )
758  {
759  converged = false;
760  xhi += 1;
761  }
762 
763  if ((xhi-xlo+1) >= hi_r)
764  {
765  converged = true;
766  }
767 
768  } while (!converged);
769 
770  /* Get precise location */
771  if (CENTERING_METHOD == CENTERING_GAUSSIAN)
772  {
773 #if WEIGHTED_FIT
774  uves_fit_1d_image(image, noise, NULL,
775 #else /* Unweighted fit like MIDAS which gives larger sigma */
776  uves_fit_1d_image(image, NULL, NULL,
777 #endif
778  true, false, false,
779  xlo, xhi, row,
780  &x0, sigma, intensity, background, slope,
781 #if WEIGHTED_FIT
782  NULL, NULL, &covariance,
783 #else
784  NULL, NULL, NULL,
785 #endif
786 
787 #if FIT_SLOPE
789 #else
791  *slope = 0;
792 #endif
793 
794  /* The fitting routine sometimes (i.e. regularly) fails
795  * because of low statistics.
796  * Recover from specific fitting errors */
797  if (cpl_error_get_code() == CPL_ERROR_NONE)
798  {
799  /* Variance is guaranteed to be positive */
800 #if WEIGHTED_FIT
801  *dx0 = sqrt(cpl_matrix_get(covariance, 0, 0));
802 #else
803  *dx0 = *sigma / sqrt(*intensity);
804 #endif
805 
806 #if WANT_BIG_LOGFILE
807  uves_msg_debug("Gaussian fit succeeded at (x, row, N) = (%f, %d, %d)",
808  x0, row, xhi-xlo+1);
809 #endif
810  }
811  else if (cpl_error_get_code() == CPL_ERROR_CONTINUE)
812  {
813  /* Fitting failed */
815 #if WANT_BIG_LOGFILE
816  uves_msg_debug("Gaussian fit failed at (x, row, N) ="
817  " (%f, %d, %d), using centroid",
818  x0, row, xhi-xlo+1);
819 #endif
820  *dx0 = *sigma / sqrt(*intensity);
821  }
822  else if (cpl_error_get_code() == CPL_ERROR_SINGULAR_MATRIX)
823  {
825 
826  /* Fitting succeeded but covariance computation failed */
827  uves_msg_debug("Covariance matrix computation failed");
828  *dx0 = *sigma / sqrt(*intensity);
829  }
830 
831  assure(cpl_error_get_code() == CPL_ERROR_NONE, cpl_error_get_code(),
832  "Gaussian fitting failed");
833 
834 #if WANT_BIG_LOGFILE
835  uves_msg_debug("Fit = (x0=%f, sigma=%f, norm=%f, backg=%f, N=%d)",
836  x0,
837  *sigma,
838  *intensity,
839  background,
840  xhi - xlo + 1);
841 #endif
842 
843  /* 'intensity' is the norm (area) of the gaussian fit.
844  But we need to return the height above zero
845  (for MIDAS compatibility).
846  height = f(x0) = background + norm/sqrt(2pi sigma^2)
847  */
848 
849  *intensity = *background + (*intensity)/(sqrt(2*M_PI) * (*sigma));
850 
851  }
852  else /* if (CENTERING_METHOD == CENTERING_GRAVITY) */
853  {
854  assure (false, CPL_ERROR_UNSUPPORTED_MODE,
855  "Centering method (no. %d) is unsupported",
856  CENTERING_METHOD);
857  }
858 
859  cleanup:
860  uves_free_matrix(&covariance);
861  return x0;
862 }
863 
int uves_gauss_linear(const double x[], const double a[], double *result)
Evaluate a gaussian with linear background.
Definition: uves_utils.c:4411
#define passure(BOOL,...)
Definition: uves_error.h:207
int uves_gauss_derivative(const double x[], const double a[], double result[])
Evaluate the derivatives of a gaussian.
Definition: uves_utils.c:4346
int uves_qclog_add_int(cpl_table *table, const char *key_name, const int value, const char *key_help, const char *format)
Add integer key to QC-LOG table.
Definition: uves_qclog.c:521
#define uves_msg(...)
Print a message on 'info' or 'debug' level.
Definition: uves_msg.h:119
int uves_gauss(const double x[], const double a[], double *result)
Evaluate a gaussian.
Definition: uves_utils.c:4291
int uves_gauss_linear_derivative(const double x[], const double a[], double result[])
Evaluate the derivatives of a gaussian with linear background.
Definition: uves_utils.c:4470
static double xcenter(const cpl_image *image, const cpl_image *noise, int xlo, int xhi, int row, centering_method CENTERING_METHOD, int bin_disp, double *sigma, double *intensity, double *dx0, double *slope, double *background)
Refine the center position of an initially detected emission line.
cpl_table * uves_wavecal_search(const cpl_image *spectrum, const cpl_image *noise, const uves_propertylist *spectrum_header, bool flat_fielded, const polynomial *order_locations, cpl_image *arcframe, int RANGE, int MINLINES, int MAXLINES, centering_method CENTERING_METHOD, int bin_disp, const int trace, const int window, cpl_table *qclog)
Search for a given number of emission lines.
double uves_polynomial_evaluate_2d(const polynomial *p, double x1, double x2)
Evaluate a 2d polynomial.
const char * uves_tostring_cpl_type(cpl_type t)
Convert a CPL type to a string.
Definition: uves_dump.c:378
#define uves_error_reset()
Definition: uves_error.h:215
#define uves_msg_debug(...)
Print a debug message.
Definition: uves_msg.h:97
#define check(CMD,...)
Definition: uves_error.h:198
double uves_pfits_get_crval2(const uves_propertylist *plist)
Find out the crval2.
Definition: uves_pfits.c:2411
static cpl_error_code detect_lines(const cpl_image *spectrum, const cpl_image *noise, const uves_propertylist *spectrum_header, bool flat_fielded, int RANGE, double THRESHOLD, centering_method CENTERING_METHOD, int bin_disp, const polynomial *order_locations, cpl_image *arcframe, cpl_table *linetable, int *ndetected, int *nrows)
Find emission lines above a certain threshold.