UVES Pipeline Reference Manual  5.4.0
uves_utils_wrappers.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: 2012-05-01 06:27:56 $
23  * $Revision: 1.123 $
24  * $Name: not supported by cvs2svn $
25  */
26 
27 #ifdef HAVE_CONFIG_H
28 # include <config.h>
29 #endif
30 
31 /*----------------------------------------------------------------------------*/
39 /*----------------------------------------------------------------------------*/
40 
41 #include <uves_utils_wrappers.h>
42 
43 #include <uves_utils.h>
44 #include <uves_utils_cpl.h>
45 #include <uves_error.h>
46 #include <uves_dump.h>
47 #include <cpl.h>
48 
49 #include <irplib_utils.h>
50 #include <stdbool.h>
51 #include <assert.h>
52 
53 /*-----------------------------------------------------------------------------
54  Local functions
55  -----------------------------------------------------------------------------*/
56 
57 static int get_candidate(const double *a, const int ia[],
58  int M, int N, int D,
59  double lambda,
60  int (*f)(const double x[], const double a[],
61  double *result),
62  int (*dfda)(const double x[], const double a[],
63  double result[]),
64  const double *x,
65  const double *y,
66  const double *sigma,
67  double *partials,
68  cpl_matrix *alpha,
69  cpl_matrix *beta,
70  double *a_da);
71 
72 static double get_chisq(int N, int D,
73  int (*f)(const double x[], const double a[],
74  double *result),
75  const double *a,
76  const double *x,
77  const double *y,
78  const double *sigma);
79 
80 
81 static cpl_image *
82 uves_image_filter_wrapper(const cpl_image *b,
83  const cpl_matrix *k,
84  cpl_filter_mode mode);
85 cpl_image *
86 uves_image_filter_median(const cpl_image * img, const cpl_matrix * mx)
87 {
88  return uves_image_filter_wrapper(img, mx, CPL_FILTER_MEDIAN);
89 }
90 
91 cpl_image *
92 uves_image_filter_linear(const cpl_image *img, const cpl_matrix * mx)
93 {
94  return uves_image_filter_wrapper(img, mx, CPL_FILTER_LINEAR);
95 
96 }
97 
98 /*----------------------------------------------------------------------------
99  Defines
100  ---------------------------------------------------------------------------*/
101 
103 #define image_is_rejected(badmap, x, y) \
104  ((badmap) != NULL && (badmap)[((x)-1) + ((y)-1)*nx] == CPL_BINARY_1)
105 
106 #ifndef UVES_FIT_MAXITER
107 #define UVES_FIT_MAXITER 1000
108 #endif
109 
110 /*-----------------------------------------------------------------------------
111  Types
112  -----------------------------------------------------------------------------*/
113 /* @cond */
114 typedef struct {
115  double x;
116  double y;
117 } uves_fit_1d_input;
118 /* @endcond */
119 
120 /*-----------------------------------------------------------------------------
121  Implementation
122  -----------------------------------------------------------------------------*/
127 static cpl_image *
128 uves_image_filter_wrapper(const cpl_image *b,
129  const cpl_matrix *k,
130  cpl_filter_mode mode)
131 {
132  const double EPSILON = 1E-5;
133  int nx = cpl_image_get_size_x(b);
134  int ny = cpl_image_get_size_y(b);
135  int nrow = cpl_matrix_get_nrow(k);
136  int ncol = cpl_matrix_get_ncol(k);
137  int i, j;
138  cpl_type type = cpl_image_get_type(b);
139  cpl_image * a = cpl_image_new(nx, ny, type);
140  // where m is a cpl_mask with a CPL_BINARY_1 whereever k has a 1.0.
141  cpl_mask* m = cpl_mask_new(ncol, nrow);
142  cpl_msg_warning(cpl_func, "nx[%d], ny[%d], ncol[%d], nrow[%d]", nx, ny, ncol, nrow);
143  for (i = 0; i < ncol ; i++)
144  {
145  for (j = 0; j < nrow ; j++)
146  {
147  double value = cpl_matrix_get(k, j, i);
148  if (fabs(value - 1.0) < EPSILON)
149  {
150  cpl_mask_set(m, i + 1, j + 1, CPL_BINARY_1);
151  }
152  }
153  }
154 
155  cpl_image_filter_mask(a, b, m, mode, CPL_BORDER_FILTER);
156  cpl_mask_delete(m);
157  return a;
158 }
159 
160 cpl_image*
161 uves_image_filter_mode(const cpl_image* b,
162  const cpl_matrix * ker,
163  cpl_filter_mode filter)
164 {
165  int nx = cpl_image_get_size_x(b);
166  int ny = cpl_image_get_size_y(b);
167  int type = cpl_image_get_type(b);
168  cpl_image * a = cpl_image_new(nx, ny, type);
169 
170  switch(filter) {
171  case CPL_FILTER_MEDIAN:
172  check_nomsg(cpl_image_filter(a, b, ker, CPL_FILTER_MEDIAN, CPL_BORDER_FILTER));
173  break;
174  case CPL_FILTER_LINEAR:
175  check_nomsg(cpl_image_filter(a, b, ker, CPL_FILTER_LINEAR, CPL_BORDER_FILTER));
176  break;
177  case CPL_FILTER_STDEV:
178  cpl_image_filter(a, b, ker, CPL_FILTER_STDEV, CPL_BORDER_FILTER);
179  break;
180  case CPL_FILTER_MORPHO:
181  cpl_image_filter(a, b, ker, CPL_FILTER_MORPHO, CPL_BORDER_FILTER);
182  break;
183  default:
184  uves_msg_error("Filter type not supported");
185  return NULL;
186  }
187  cleanup:
188 
189  return a;
190 
191 }
192 
193 
194 /*----------------------------------------------------------------------------*/
201 /*----------------------------------------------------------------------------*/
202 void uves_image_reject_all(cpl_image *image)
203 {
204  int nx, ny;
205  int x, y;
206 
207  assure_nomsg( image != NULL, CPL_ERROR_NULL_INPUT );
208  nx = cpl_image_get_size_x(image);
209  ny = cpl_image_get_size_y(image);
210 
211  for (y = 1; y <= ny; y++) {
212  for (x = 1; x <= nx; x++) {
213  cpl_image_reject(image, x, y);
214  }
215  }
216 
217  cleanup:
218  return;
219 }
220 
221 
222 /*----------------------------------------------------------------------------*/
281 /*----------------------------------------------------------------------------*/
282 
283 cpl_error_code
284 uves_fit_1d_image(const cpl_image *image, const cpl_image *noise,
285  const cpl_binary *image_badmap,
286  bool horizontal, bool fix_back, bool fit_back,
287  int xlo, int xhi, int y_0,
288  double *x0, double *sigma, double *norm, double *background,
289  double *slope,
290  double *mse, double *red_chisq,
291  cpl_matrix **covariance,
292  int (*f) (const double x[], const double a[], double *result),
293  int (*dfda)(const double x[], const double a[], double result[]),
294  int M)
295 {
296  cpl_vector *x = NULL;
297  cpl_vector *y = NULL;
298  cpl_vector *sigma_y = NULL; /* Noise vector */
299 
300  cpl_fit_mode fit_pattern;
301  int nx, ny; /* Image size */
302  int N; /* Number of good pixels */
303  int i, j;
304  cpl_type image_type;
305 
306  const double *image_data = NULL; /* Pointer to data */
307  const double *noise_data = NULL; /* Pointer to data */
308 
309  assure( x0 != NULL , CPL_ERROR_NULL_INPUT, "Null fit parameter");
310  assure( sigma != NULL , CPL_ERROR_NULL_INPUT, "Null fit parameter");
311  assure( norm != NULL , CPL_ERROR_NULL_INPUT, "Null fit parameter");
312  assure( background != NULL, CPL_ERROR_NULL_INPUT, "Null fit parameter");
313  /* mse, red_chisq, covariance may be NULL */
314 
315  assure( image != NULL, CPL_ERROR_NULL_INPUT, "Null image");
316 
317  image_type = cpl_image_get_type(image);
318 
319  /* To support the following types, use cpl_image_get() or
320  more multiple type pointers to data buffer.
321  cpl_ensure_code( image_type == CPL_TYPE_INT ||
322  image_type == CPL_TYPE_FLOAT ||
323  image_type == CPL_TYPE_DOUBLE,
324  CPL_ERROR_TYPE_MISMATCH);
325  */
326 
327  assure( image_type == CPL_TYPE_DOUBLE, CPL_ERROR_UNSUPPORTED_MODE,
328  "Unsupported type: %s", uves_tostring_cpl_type(image_type));
329 
330  image_data = cpl_image_get_data_double_const(image);
331 
332  if (noise != NULL)
333  {
334  image_type = cpl_image_get_type(noise);
335  /* See comment above.
336  cpl_ensure_code( image_type == CPL_TYPE_INT ||
337  image_type == CPL_TYPE_FLOAT ||
338  image_type == CPL_TYPE_DOUBLE,
339  CPL_ERROR_TYPE_MISMATCH);
340  */
341  assure( image_type == CPL_TYPE_DOUBLE, CPL_ERROR_UNSUPPORTED_MODE,
342  "Unsupported type: %s", uves_tostring_cpl_type(image_type));
343 
344  noise_data = cpl_image_get_data_double_const(noise);
345  }
346 
347  nx = cpl_image_get_size_x(image);
348  ny = cpl_image_get_size_y(image);
349 
350  if (horizontal)
351  {
352  assure( 1 <= xlo && xlo <= xhi && xhi <= nx &&
353  1 <= y_0 && y_0 <= ny, CPL_ERROR_ACCESS_OUT_OF_RANGE,
354  "Illegal window (%d, %d)-(%d, %d), image: %dx%d",
355  xlo, y_0, xhi, y_0,
356  nx, ny);
357  }
358  else
359  {
360  assure( 1 <= xlo && xlo <= xhi && xhi <= ny &&
361  1 <= y_0 && y_0 <= nx,
362  CPL_ERROR_ACCESS_OUT_OF_RANGE,
363  "Illegal window (%d, %d)-(%d, %d), image: %dx%d",
364  y_0, xlo, y_0, xhi,
365  nx, ny);
366  }
367 
368  /* Noise image must have same size
369  * as the input image
370  */
371  if (noise != NULL)
372  {
373  assure( cpl_image_get_size_x(noise) == nx &&
374  cpl_image_get_size_y(noise) == ny,
375  CPL_ERROR_INCOMPATIBLE_INPUT, "Noise image: %" CPL_SIZE_FORMAT "x%" CPL_SIZE_FORMAT ", image: %dx%d:",
376  cpl_image_get_size_x(noise),
377  cpl_image_get_size_y(noise),
378  nx, ny);
379  }
380 
381  /* Count good pixels in sub-window, check that noise image is positive */
382  N = 0;
383  for (i = xlo; i <= xhi; i++)
384  {
385  if ( !image_is_rejected(image_badmap,
386  (horizontal) ? i : y_0,
387  (horizontal) ? y_0 : i) )
388  {
389  if ( noise != NULL)
390  {
391  if( !image_is_rejected(image_badmap,
392  (horizontal) ? i : y_0,
393  (horizontal) ? y_0 : i) )
394  {
395  /* Noise image must be positive (only check
396  pixels that are actually used) */
397  assure( /*cpl_image_get(noise,
398  (horizontal) ? i : y_0,
399  (horizontal) ? y_0 : i,
400  &pis_rejected)*/
401  noise_data[(((horizontal) ? i : y_0) - 1) +
402  (((horizontal) ? y_0 : i) - 1) * nx]
403  > 0,
404  CPL_ERROR_ILLEGAL_INPUT,
405  "Non-positive noise at (%d, %d): %e",
406  (horizontal) ? i : y_0,
407  (horizontal) ? y_0 : i,
408  noise_data[(((horizontal) ? i : y_0) - 1) +
409  (((horizontal) ? y_0 : i) - 1) * nx]
410  );
411 
412  N += 1;
413  }
414  else
415  {
416  /* Pixel value is good, but noise pixel is
417  bad. Don't use. */
418  }
419  }
420  else
421  {
422  /* Pixel is good. No noise image */
423  N += 1;
424  }
425  }
426  }
427 
428  /* Check that there is at least one good pixel */
429  assure( N >= 1, CPL_ERROR_ILLEGAL_INPUT, "Only %d good pixel(s)", N);
430 
431  /* Allocate space */
432  x = cpl_vector_new(N);
433  y = cpl_vector_new(N);
434  if (noise != NULL)
435  {
436  sigma_y = cpl_vector_new(N);
437  assure_mem( sigma_y );
438  }
439 
440  if (fix_back)
441  {
442  fit_pattern = CPL_FIT_CENTROID | CPL_FIT_STDEV | CPL_FIT_AREA;
443  }
444  else if (fit_back)
445  {
446  fit_pattern = CPL_FIT_AREA | CPL_FIT_OFFSET;
447  }
448  else
449  {
450  fit_pattern = CPL_FIT_ALL;
451  }
452 
453  assure_mem( x );
454  assure_mem( y );
455 
456  /* Copy the N good pixels from the input image to vectors,
457  j count good pixels */
458  for (i = xlo, j = 0;
459  i <= xhi;
460  i++)
461  {
462  double flux;
463 
464  /*
465  flux = cpl_image_get(image,
466  (horizontal) ? xlo+i : y_0,
467  (horizontal) ? y_0 : xlo+i,
468  &pis_rejected);
469  */
470 
471  flux = image_data[(((horizontal) ? i : y_0) - 1) +
472  (((horizontal) ? y_0 : i) - 1) * nx];
473 
474  //if (!pis_rejected)
475  if ( !image_is_rejected(image_badmap,
476  (horizontal) ? i : y_0,
477  (horizontal) ? y_0 : i) )
478  {
479  if (noise != NULL)
480  {
481  double flux_noise;
482 
483  /* flux_noise = cpl_image_get(noise,
484  (horizontal) ? xlo+i : y_0,
485  (horizontal) ? y_0 : xlo+i,
486  &pis_rejected);
487  */
488 
489  flux_noise = noise_data
490  [(((horizontal) ? i : y_0) - 1) +
491  (((horizontal) ? y_0 : i) - 1) * nx];
492 
493  //if (!pis_rejected)
494  if ( !image_is_rejected(image_badmap,
495  (horizontal) ?
496  i : y_0,
497  (horizontal)
498  ? y_0 : i) )
499  {
500  cpl_vector_set(x, j, i);
501  cpl_vector_set(y, j, flux);
502  cpl_vector_set(sigma_y, j, flux_noise);
503  j++;
504  }
505  }
506  else
507  {
508  cpl_vector_set(x, j, i);
509  cpl_vector_set(y, j, flux);
510  j++;
511  }
512  }
513  }
514  passure( j == N, "%d %d", j, N);
515 
516  check( uves_fit_1d(x, NULL, /* x, sigma_x */
517  y, sigma_y,
518  fit_pattern, fit_back,
519  x0, sigma, norm, background,
520  slope,
521  mse, red_chisq,
522  covariance,
523  f, dfda, M),
524  "Fit failed");
525 
526 
527  cleanup:
528  uves_free_vector(&x);
529  uves_free_vector(&y);
530  uves_free_vector(&sigma_y);
531 
532  return cpl_error_get_code();
533 }
534 
535 /*----------------------------------------------------------------------------*/
543 /*----------------------------------------------------------------------------*/
544 static int uves_fit_1d_compare(const void *left,
545  const void *right)
546 {
547  return
548  (((uves_fit_1d_input *)left )->x <
549  ((uves_fit_1d_input *)right)->x) ? -1 :
550  (((uves_fit_1d_input *)left )->x ==
551  ((uves_fit_1d_input *)right)->x) ? 0 : 1;
552 }
553 
554 /*----------------------------------------------------------------------------*/
577 /*----------------------------------------------------------------------------*/
578 cpl_error_code
579 uves_fit_1d(cpl_vector *x, const cpl_vector *sigma_x,
580  cpl_vector *y, const cpl_vector *sigma_y,
581  cpl_fit_mode fit_pars, bool fit_back,
582  double *x0, double *sigma, double *area, double *offset, double *slope,
583  double *mse, double *red_chisq,
584  cpl_matrix **covariance,
585  int (*f) (const double x[], const double a[], double *result),
586  int (*dfda)(const double x[], const double a[], double result[]),
587  int M)
588 {
589  cpl_matrix *x_matrix = NULL; /* LM algorithm needs a matrix,
590  not a vector */
591 
592  int N; /* Number of data points */
593  double xlo, xhi; /* Min/max x */
594 
595  /* Initial parameter values */
596  double x0_guess = 0; /* Avoid warnings about uninitialized variables */
597  double sigma_guess = 0;
598  double area_guess;
599  double offset_guess;
600 
601  cpl_ensure_code( M == 4 || M == 5, CPL_ERROR_UNSUPPORTED_MODE);
602 
603  /* Validate input */
604  cpl_ensure_code( x != NULL, CPL_ERROR_NULL_INPUT);
605  cpl_ensure_code( sigma_x == NULL, CPL_ERROR_UNSUPPORTED_MODE);
606  cpl_ensure_code( y != NULL, CPL_ERROR_NULL_INPUT);
607  /* sigma_y may be NULL or non-NULL */
608 
609  N = cpl_vector_get_size(x);
610 
611  cpl_ensure_code( N == cpl_vector_get_size(y),
612  CPL_ERROR_INCOMPATIBLE_INPUT);
613 
614  if (sigma_x != NULL)
615  {
616  cpl_ensure_code( N == cpl_vector_get_size(sigma_x),
617  CPL_ERROR_INCOMPATIBLE_INPUT);
618  }
619  if (sigma_y != NULL)
620  {
621  cpl_ensure_code( N == cpl_vector_get_size(sigma_y),
622  CPL_ERROR_INCOMPATIBLE_INPUT);
623  }
624 
625  cpl_ensure_code( x0 != NULL &&
626  sigma != NULL &&
627  area != NULL &&
628  offset != NULL &&
629  (M != 5 || slope != NULL), CPL_ERROR_NULL_INPUT);
630 
631  if (! (fit_pars & CPL_FIT_STDEV))
632  {
633  cpl_ensure_code( *sigma > 0, CPL_ERROR_ILLEGAL_INPUT);
634  }
635 
636  cpl_ensure_code( !fit_back || fit_pars == (CPL_FIT_OFFSET | CPL_FIT_AREA),
637  CPL_ERROR_INCOMPATIBLE_INPUT);
638 
639  /* Input area must be positive if fit_back is false */
640  if (! (fit_pars & CPL_FIT_AREA) && !fit_back)
641  {
642  cpl_ensure_code( *area > 0, CPL_ERROR_ILLEGAL_INPUT);
643  }
644 
645  /* mse, red_chisq may be NULL */
646 
647  /* Need more than number_of_parameters points to calculate chi^2.
648  * There are less than 5 parameters. */
649  cpl_ensure_code( red_chisq == NULL || N >= 5, CPL_ERROR_ILLEGAL_INPUT);
650 
651  if (covariance != NULL) *covariance = NULL;
652  /* If covariance computation is requested, then
653  * return either the covariance matrix or NULL
654  * (don't return undefined pointer).
655  */
656 
657  /* Cannot compute chi square & covariance without sigma_y */
658  cpl_ensure_code( (red_chisq == NULL && covariance == NULL) ||
659  sigma_y != NULL,
660  CPL_ERROR_INCOMPATIBLE_INPUT);
661 
662  /* Create matrix from x-data */
663  x_matrix = cpl_matrix_wrap(N, 1, cpl_vector_get_data(x));
664  if (x_matrix == NULL)
665  {
666  cpl_ensure_code(
667  CPL_FALSE,
668  CPL_ERROR_ILLEGAL_OUTPUT);
669  }
670 
671  /* Check that any provided sigmas are positive. */
672  if (sigma_x != NULL && cpl_vector_get_min(sigma_x) <= 0)
673  {
674  cpl_matrix_unwrap(x_matrix);
675  cpl_ensure_code(
676  CPL_FALSE,
677  CPL_ERROR_ILLEGAL_INPUT);
678  }
679  if (sigma_y != NULL && cpl_vector_get_min(sigma_y) <= 0)
680  {
681  cpl_matrix_unwrap(x_matrix);
682  cpl_ensure_code(
683  CPL_FALSE,
684  CPL_ERROR_ILLEGAL_INPUT);
685  }
686 
687  /* Compute guess parameters using robust estimation
688  * (This step might be improved by taking into account the
689  * uncertainties but for simplicity's sake do not)
690  */
691  if (fit_back)
692  {
693  /* We need to estimate only these two parameters */
694  assert( fit_pars == CPL_FIT_OFFSET || CPL_FIT_AREA);
695 
696  offset_guess = cpl_vector_get_median_const(y);
697  area_guess = N*(cpl_vector_get_mean(y) - offset_guess);
698  /* Sum of (flux-offset) */
699 
700  x0_guess = *x0;
701  sigma_guess = *sigma;
702  }
703  else {
704  double sum = 0.0;
705  double quartile[3];
706  double fraction[3] = {0.25 , 0.50 , 0.75};
707  const double *y_data = cpl_vector_get_data_const(y);
708 
709  if (fit_pars & CPL_FIT_OFFSET)
710  {
711  /* Estimate offset as 25% percentile of y-values.
712  (The minimum value may be too low for noisy input,
713  the median is too high if there is not much
714  background in the supplied data, so use
715  something inbetween).
716  */
717 
718  cpl_vector *y_dup = cpl_vector_duplicate(y);
719 
720  if (y_dup == NULL)
721  {
722  cpl_matrix_unwrap(x_matrix);
723  cpl_ensure_code(
724  CPL_FALSE,
725  CPL_ERROR_ILLEGAL_OUTPUT);
726  }
727 
728  offset_guess = uves_utils_get_kth_double(
729  cpl_vector_get_data(y_dup), N, N/4);
730 
731  cpl_vector_delete(y_dup);
732  }
733  else
734  {
735  offset_guess = *offset;
736  }
737 
738  /* Get quartiles of distribution
739  (only bother if it's needed for estimation of x0 or sigma) */
740  if ( (fit_pars & CPL_FIT_CENTROID) ||
741  (fit_pars & CPL_FIT_STDEV )
742  )
743  {
744  /* The algorithm requires the input to be sorted
745  as function of x, so do that (using qsort), and
746  work on the sorted copy. Of course, the y-vector
747  must be re-ordered along with x.
748  sigma_x and sigma_y are not used, so don't copy those.
749  */
750 
751  uves_fit_1d_input
752  *sorted_input = cpl_malloc(N * sizeof(*sorted_input));
753  const double *x_data = cpl_matrix_get_data_const(x_matrix);
754  cpl_boolean is_sorted = CPL_TRUE;
755  int i;
756 
757  if (sorted_input == NULL)
758  {
759  cpl_matrix_unwrap(x_matrix);
760  cpl_ensure_code(
761  CPL_FALSE,
762  CPL_ERROR_ILLEGAL_INPUT);
763  }
764 
765  for (i = 0; i < N; i++)
766  {
767  sorted_input[i].x = x_data[i];
768  sorted_input[i].y = y_data[i];
769 
770  is_sorted = is_sorted &&
771  (i==0 || (x_data[i-1] < x_data[i]));
772  }
773 
774  if (!is_sorted)
775  {
776  qsort(sorted_input, N, sizeof(*sorted_input),
777  &uves_fit_1d_compare);
778  }
779 
780  for(i = 0; i < N; i++)
781  {
782  double flux = sorted_input[i].y;
783 
784  sum += (flux - offset_guess);
785  }
786  /* Note that 'sum' must be calculated the same way as
787  'running_sum' below, Otherwise (due to round-off error)
788  'running_sum' might end up being different from 'sum'(!).
789  Specifically, it will not work to calculate 'sum' as
790 
791  (flux1 + ... + fluxN) - N*offset_guess
792  */
793 
794  if (sum > 0.0)
795  {
796  double flux, x1, x2;
797  double running_sum = 0.0;
798  int j;
799 
800  i = 0;
801  flux = sorted_input[i].y - offset_guess;
802 
803  for (j = 0; j < 3; j++)
804  {
805  double limit = fraction[j] * sum;
806  double k;
807 
808  while (running_sum + flux < limit && i < N-1)
809  {
810  running_sum += flux;
811  i++;
812  flux = sorted_input[i].y - offset_guess;
813  }
814 
815  /* Fraction [0;1] of current flux needed
816  to reach the quartile */
817  k = (limit - running_sum)/flux;
818 
819  if (k <= 0.5)
820  {
821  /* Interpolate linearly between
822  current and previous position
823  */
824  if (0 < i)
825  {
826  x1 = sorted_input[i-1].x;
827  x2 = sorted_input[i].x;
828 
829  quartile[j] =
830  x1*(0.5-k) +
831  x2*(0.5+k);
832  /*
833  k=0 => quartile = midpoint,
834  k=0.5 => quartile = x2
835  */
836  }
837  else
838  {
839  quartile[j] = sorted_input[i].x;
840  }
841  }
842  else
843  {
844  /* Interpolate linearly between
845  current and next position */
846  if (i < N-1)
847  {
848  x1 = sorted_input[i].x;
849  x2 = sorted_input[i+1].x;
850 
851  quartile[j] =
852  x1*( 1.5-k) +
853  x2*(-0.5+k);
854  /*
855  k=0.5 => quartile = x1,
856  k=1.0 => quartile = midpoint
857  */
858  }
859  else
860  {
861  quartile[j] = sorted_input[i].x;
862  }
863  }
864  }
865  }
866  else
867  {
868  /* If there's no flux (sum = 0) then
869  set quartiles to something that's not
870  completely insensible.
871  */
872  quartile[1] = cpl_matrix_get_mean(x_matrix);
873 
874  quartile[2] = quartile[1];
875  quartile[0] = quartile[2] - 1.0;
876  /* Then sigma_guess = 1.0 */
877  }
878 
879  cpl_free(sorted_input);
880  } /* If need to compute quartiles */
881 
882  /* x0_guess = median of distribution */
883  x0_guess = (fit_pars & CPL_FIT_CENTROID) ? quartile[1] : *x0;
884 
885  /* sigma_guess = median of absolute residuals
886  *
887  * (68% is within 1 sigma, and 50% is within 0.6744
888  * sigma, => quartile3-quartile1 = 2*0.6744 sigma)
889  */
890  sigma_guess = (fit_pars & CPL_FIT_STDEV) ?
891  (quartile[2] - quartile[0]) / (2*0.6744) : *sigma;
892 
893  area_guess = (fit_pars & CPL_FIT_AREA) ?
894  (cpl_vector_get_max(y) - offset_guess) * sqrt(2*M_PI) * sigma_guess : *area;
895  /* This formula makes sense only if the area is positive */
896 
897  /* Make sure that the area is a positive number */
898  if ( area_guess <= 0) area_guess = 1.0;
899  if (sigma_guess <= 0) sigma_guess = 1.0;
900  }
901 
902  /* Wrap parameters, fit, unwrap */
903  {
904  cpl_vector *a;
905  int ia[5]; /* The last element
906  is ignored if
907  M = 4 */
908  cpl_error_code ec;
909 
910  assert(M == 4 || M == 5);
911  a = cpl_vector_new(M);
912 
913  if (a == NULL)
914  {
915  cpl_matrix_unwrap(x_matrix);
916  cpl_ensure_code(
917  CPL_FALSE,
918  CPL_ERROR_ILLEGAL_OUTPUT);
919  }
920 
921  cpl_vector_set(a, 0, x0_guess);
922  cpl_vector_set(a, 1, sigma_guess);
923  cpl_vector_set(a, 2, area_guess);
924  cpl_vector_set(a, 3, offset_guess);
925 
926  ia[0] = fit_pars & CPL_FIT_CENTROID;
927  ia[1] = fit_pars & CPL_FIT_STDEV;
928  ia[2] = fit_pars & CPL_FIT_AREA;
929  ia[3] = fit_pars & CPL_FIT_OFFSET;
930 
931  if (M == 5)
932  {
933  /* linear sky-term, first hold it constant,
934  * then call LM-fitting a second time where
935  * it is non-constant */
936  if (fit_pars & CPL_FIT_OFFSET)
937  {
938  cpl_vector_set(a, 4, 0);
939  }
940  else
941  {
942  cpl_vector_set(a, 4, *slope);
943  }
944 
945  ia[4] = 0;
946  }
947 
948  ec = uves_fit(x_matrix, NULL,
949  y, sigma_y,
950  a, ia, f, dfda,
951  mse, red_chisq,
952  covariance);
953 
954 // printf("Sky: e='%s'\n", cpl_error_get_message()); cpl_vector_dump(a, stdout);
955  if (M == 5)
956  {
957  ia[4] = fit_pars & CPL_FIT_OFFSET;
958 
959  if (covariance != NULL)
960  {
961  uves_free_matrix(covariance);
962  }
963 
964  ec = uves_fit(x_matrix, NULL,
965  y, sigma_y,
966  a, ia, f, dfda,
967  mse, red_chisq,
968  covariance);
969 // printf("Sky: e='%s'\n", cpl_error_get_message()); cpl_vector_dump(a, stdout);
970  }
971 
972  cpl_matrix_unwrap(x_matrix);
973 
974  /* Check return status of fitting routine */
975  if (ec == CPL_ERROR_NONE ||
976  ec == CPL_ERROR_SINGULAR_MATRIX)
977  {
978  /* The LM algorithm converged. The computation
979  * of the covariance matrix might have failed.
980  */
981 
982  /* In principle, the LM algorithm might have converged
983  * to a negative sigma (even if the guess value was
984  * positive). Make sure that the returned sigma is positive
985  * (by convention).
986  */
987 
988  if (CPL_FIT_CENTROID) *x0 = cpl_vector_get(a, 0);
989  if (CPL_FIT_STDEV ) *sigma = fabs(cpl_vector_get(a, 1));
990  if (CPL_FIT_AREA ) *area = cpl_vector_get(a, 2);
991  if (CPL_FIT_OFFSET ) {
992  *offset = cpl_vector_get(a, 3);
993  if (M == 5) *slope = cpl_vector_get(a, 4);
994  }
995  }
996 
997  cpl_vector_delete(a);
998 
999  xlo = cpl_vector_get_min(x);
1000  xhi = cpl_vector_get_max(x);
1001 
1002  if (ec == CPL_ERROR_CONTINUE ||
1003  !(
1004  !irplib_isnan(*x0 ) && !irplib_isinf(*x0 ) &&
1005  !irplib_isnan(*sigma ) && !irplib_isinf(*sigma ) &&
1006  !irplib_isnan(*area ) && !irplib_isinf(*area ) &&
1007  !irplib_isnan(*offset) && !irplib_isinf(*offset) &&
1008  ((M != 5) || (!irplib_isnan(*slope ) && !irplib_isinf(*slope ))) &&
1009  xlo <= *x0 && *x0 <= xhi &&
1010  0 < *sigma && *sigma < (xhi - xlo + 1) &&
1011  (fit_back || 0 < *area)
1012  /* This extra check on the background level makes sense
1013  iff the input flux is assumed to be positive
1014  && *offset > - *area */
1015  )
1016  )
1017  {
1018  /* The LM algorithm did not converge, or it converged to
1019  * a non-sensical result. Return the guess parameter values
1020  * in order to enable the caller to recover. */
1021 
1022  *x0 = x0_guess;
1023  *sigma = sigma_guess;
1024  *area = area_guess;
1025  *offset = offset_guess;
1026  if (M == 5) *slope = 0;
1027 
1028  /* In this case the covariance matrix will not make sense
1029  (because the LM algorithm failed), so delete it */
1030  if (covariance != NULL && *covariance != NULL)
1031  {
1032  cpl_matrix_delete(*covariance);
1033  *covariance = NULL;
1034  }
1035 
1036  /* Return CPL_ERROR_CONTINUE if the fitting routine failed */
1037  cpl_ensure_code(
1038  CPL_FALSE,
1039  CPL_ERROR_CONTINUE);
1040  }
1041  }
1042 
1043  return CPL_ERROR_NONE;
1044 }
1045 /* @endcond */
1046 
1047 #define DEBUG_LM 0 /* Set to non-zero to print info on the error msg level */
1048 /*----------------------------------------------------------------------------*/
1115 /*----------------------------------------------------------------------------*/
1116 cpl_error_code
1117 uves_fit(const cpl_matrix *x, const cpl_matrix *sigma_x,
1118  const cpl_vector *y, const cpl_vector *sigma_y,
1119  cpl_vector *a, const int ia[],
1120  int (*f)(const double x[], const double a[], double *result),
1121  int (*dfda)(const double x[], const double a[], double result[]),
1122  double *mse,
1123  double *red_chisq,
1124  cpl_matrix **covariance)
1125 {
1126  const double *x_data = NULL; /* Pointer to input data */
1127  const double *y_data = NULL; /* Pointer to input data */
1128  const double *sigma_data = NULL; /* Pointer to input data */
1129  int N = 0; /* Number of data points */
1130  int D = 0; /* Dimension of x-points */
1131  int M = 0; /* Number of fit parameters */
1132  int Mfit = 0; /* Number of non-constant fit
1133  parameters */
1134 
1135  double lambda = 0.0; /* Lambda in L-M algorithm */
1136  double MAXLAMBDA = 10e40; /* Parameter to control the graceful exit
1137  if steepest descent unexpectedly fails */
1138  double chi_sq = 0.0; /* Current chi^2 */
1139  int count = 0; /* Number of successive small improvements
1140  in chi^2 */
1141  int iterations = 0;
1142 
1143  cpl_matrix *alpha = NULL; /* The MxM ~curvature matrix used in L-M */
1144  cpl_matrix *beta = NULL; /* Mx1 matrix = -.5 grad(chi^2) */
1145  double *a_data = NULL; /* Parameters, a */
1146  double *a_da = NULL; /* Candidate position a+da */
1147  double *part = NULL; /* The partial derivatives df/da */
1148  int *ia_local = NULL; /* non-NULL version of ia */
1149 
1150  /* If covariance computation is requested, then either
1151  * return the covariance matrix or return NULL.
1152  */
1153  if (covariance != NULL) *covariance = NULL;
1154 
1155  /* Validate input */
1156  cpl_ensure_code(x != NULL, CPL_ERROR_NULL_INPUT);
1157  cpl_ensure_code(sigma_x == NULL, CPL_ERROR_UNSUPPORTED_MODE);
1158  cpl_ensure_code(y != NULL, CPL_ERROR_NULL_INPUT);
1159  cpl_ensure_code(a != NULL, CPL_ERROR_NULL_INPUT);
1160  /* ia may be NULL */
1161  cpl_ensure_code(f != NULL, CPL_ERROR_NULL_INPUT);
1162  cpl_ensure_code(dfda != NULL, CPL_ERROR_NULL_INPUT);
1163 
1164  /* Chi^2 and covariance computations require sigmas to be known */
1165  cpl_ensure_code( sigma_y != NULL || (red_chisq == NULL && covariance == NULL),
1166  CPL_ERROR_INCOMPATIBLE_INPUT);
1167 
1168  D = cpl_matrix_get_ncol(x);
1169  N = cpl_matrix_get_nrow(x);
1170  M = cpl_vector_get_size(a);
1171  cpl_ensure_code(N > 0 && D > 0 && M > 0, CPL_ERROR_ILLEGAL_INPUT);
1172 
1173  cpl_ensure_code( cpl_vector_get_size(y) == N,
1174  CPL_ERROR_INCOMPATIBLE_INPUT);
1175 
1176  x_data = cpl_matrix_get_data_const(x);
1177  y_data = cpl_vector_get_data_const(y);
1178  a_data = cpl_vector_get_data(a);
1179 
1180  if (sigma_y != NULL)
1181  {
1182  cpl_ensure_code( cpl_vector_get_size(sigma_y) == N,
1183  CPL_ERROR_INCOMPATIBLE_INPUT);
1184  /* Sigmas must be positive */
1185  cpl_ensure_code( cpl_vector_get_min (sigma_y) > 0,
1186  CPL_ERROR_ILLEGAL_INPUT);
1187  sigma_data = cpl_vector_get_data_const(sigma_y);
1188  }
1189 
1190  ia_local = cpl_malloc(M * sizeof(int));
1191  cpl_ensure_code(ia_local != NULL, CPL_ERROR_ILLEGAL_OUTPUT);
1192 
1193  /* Count non-constant fit parameters, copy ia */
1194  if (ia != NULL)
1195  {
1196  int i;
1197 
1198  Mfit = 0;
1199  for (i = 0; i < M; i++)
1200  {
1201  ia_local[i] = ia[i];
1202 
1203  if (ia[i] != 0)
1204  {
1205  Mfit += 1;
1206  }
1207  }
1208 
1209  if (! (Mfit > 0))
1210  {
1211  cpl_free(ia_local);
1212  cpl_ensure_code( CPL_FALSE,
1213  CPL_ERROR_ILLEGAL_INPUT);
1214  }
1215  }
1216  else
1217  {
1218  /* All parameters participate */
1219  int i;
1220 
1221  Mfit = M;
1222 
1223  for (i = 0; i < M; i++)
1224  {
1225  ia_local[i] = 1;
1226  }
1227  }
1228 
1229  /* To compute reduced chi^2, we need N > Mfit */
1230  if (! ( red_chisq == NULL || N > Mfit ) )
1231  {
1232  cpl_free(ia_local);
1233  cpl_ensure_code(
1234  CPL_FALSE,
1235  CPL_ERROR_ILLEGAL_INPUT);
1236  }
1237 
1238  /* Create alpha, beta, a_da, part work space */
1239  alpha = cpl_matrix_new(Mfit, Mfit);
1240  if (alpha == NULL)
1241  {
1242  cpl_free(ia_local);
1243  cpl_ensure_code(
1244  CPL_FALSE,
1245  CPL_ERROR_ILLEGAL_OUTPUT);
1246  }
1247 
1248  beta = cpl_matrix_new(Mfit, 1);
1249  if (beta == NULL)
1250  {
1251  cpl_free(ia_local);
1252  cpl_matrix_delete(alpha);
1253  cpl_ensure_code(
1254  CPL_FALSE,
1255  CPL_ERROR_ILLEGAL_OUTPUT);
1256  }
1257 
1258  a_da = cpl_malloc(M * sizeof(double));
1259  if (a_da == NULL)
1260  {
1261  cpl_free(ia_local);
1262  cpl_matrix_delete(alpha);
1263  cpl_matrix_delete(beta);
1264  cpl_ensure_code(
1265  CPL_FALSE,
1266  CPL_ERROR_ILLEGAL_OUTPUT);
1267  }
1268 
1269  part = cpl_malloc(M * sizeof(double));
1270  if (part == NULL)
1271  {
1272  cpl_free(ia_local);
1273  cpl_matrix_delete(alpha);
1274  cpl_matrix_delete(beta);
1275  cpl_free(a_da);
1276  cpl_ensure_code(
1277  CPL_FALSE,
1278  CPL_ERROR_ILLEGAL_OUTPUT);
1279  }
1280 
1281  /* Initialize loop variables */
1282  lambda = 0.001;
1283  count = 0;
1284  iterations = 0;
1285  if( (chi_sq = get_chisq(N, D, f, a_data, x_data, y_data, sigma_data)) < 0)
1286  {
1287  cpl_free(ia_local);
1288  cpl_matrix_delete(alpha);
1289  cpl_matrix_delete(beta);
1290  cpl_free(a_da);
1291  cpl_free(part);
1292  cpl_ensure_code(
1293  CPL_FALSE,
1294  cpl_error_get_code());
1295  }
1296 
1297 #if DEBUG_LM
1298  uves_msg_error("Initial chi^2 = %f", chi_sq);
1299  {int i;
1300  for (i = 0; i < M; i++)
1301  {
1302  uves_msg_error("Initial a[%d] = %e", i, a_data[i]);
1303  }
1304  }
1305 #endif
1306 
1307  /* Iterate until chi^2 didn't improve substantially many (say, 5)
1308  times in a row */
1309  while (count < 5 && lambda < MAXLAMBDA && iterations < UVES_FIT_MAXITER)
1310  {
1311  /* In each iteration lambda increases, or chi^2 decreases or
1312  count increases. Because chi^2 is bounded from below
1313  (and lambda and count from above), the loop will terminate */
1314 
1315  double chi_sq_candidate = 0.0;
1316  int returncode = 0;
1317 
1318  /* Get candidate position in parameter space = a+da,
1319  * where alpha * da = beta .
1320  * Increase lambda until alpha is non-singular
1321  */
1322 
1323  while( (returncode = get_candidate(a_data, ia_local,
1324  M, N, D,
1325  lambda, f, dfda,
1326  x_data, y_data, sigma_data,
1327  part, alpha, beta, a_da)
1328  ) != 0
1329  && cpl_error_get_code() == CPL_ERROR_SINGULAR_MATRIX
1330  && lambda < MAXLAMBDA)
1331  {
1332 #if DEBUG_LM
1333  uves_msg_error("Singular matrix. lambda = %e", lambda);
1334 #endif
1335  cpl_error_reset();
1336  lambda *= 9.0;
1337  }
1338 
1339  /* Set error if lambda diverged */
1340  if ( !( lambda < MAXLAMBDA ) )
1341  {
1342  cpl_free(ia_local);
1343  cpl_matrix_delete(alpha);
1344  cpl_matrix_delete(beta);
1345  cpl_free(a_da);
1346  cpl_free(part);
1347  cpl_ensure_code(
1348  CPL_FALSE,
1349  CPL_ERROR_CONTINUE);
1350  }
1351 
1352  if (returncode != 0)
1353  {
1354  cpl_free(ia_local);
1355  cpl_matrix_delete(alpha);
1356  cpl_matrix_delete(beta);
1357  cpl_free(a_da);
1358  cpl_free(part);
1359  cpl_ensure_code(
1360  CPL_FALSE,
1361  cpl_error_get_code());
1362  }
1363 
1364  /* Get chi^2(a+da) */
1365  if ((chi_sq_candidate = get_chisq(N, D, f, a_da,
1366  x_data, y_data, sigma_data)) < 0)
1367  {
1368  cpl_free(ia_local);
1369  cpl_matrix_delete(alpha);
1370  cpl_matrix_delete(beta);
1371  cpl_free(a_da);
1372  cpl_free(part);
1373  cpl_ensure_code(
1374  CPL_FALSE,
1375  cpl_error_get_code());
1376  }
1377 
1378  if (chi_sq_candidate > 1.000001 * chi_sq)
1379  {
1380  /* Move towards steepest descent */
1381 #if DEBUG_LM
1382  uves_msg_error("Chi^2 = %f Candidate = %f Lambda = %e",
1383  chi_sq, chi_sq_candidate, lambda);
1384 #endif
1385  lambda *= 9.0;
1386  }
1387  else
1388  {
1389 #if DEBUG_LM
1390  uves_msg_error("Chi^2 = %f Candidate = %f* Lambda = %e count = %d",
1391  chi_sq, chi_sq_candidate, lambda, count);
1392 #endif
1393 
1394  /* Move towards Newton's algorithm */
1395  lambda /= 10.0;
1396 
1397  /* Count the number of successive improvements in chi^2 of
1398  less than 0.01 relative */
1399  if ( chi_sq == 0 ||
1400  (chi_sq - chi_sq_candidate)/chi_sq < .01)
1401  {
1402  count += 1;
1403  }
1404  else
1405  {
1406  /* Chi^2 improved by a significant amount,
1407  reset counter */
1408  count = 0;
1409  }
1410 
1411  if (chi_sq_candidate < chi_sq)
1412  {
1413  /* chi^2 improved, update a and chi^2 */
1414  int i;
1415  for (i = 0; i < M; i++)
1416  {
1417  a_data[i] = a_da[i];
1418 #if DEBUG_LM
1419  uves_msg_error("-> a[%d] = %e", i, a_da[i]);
1420 #endif
1421  }
1422  chi_sq = chi_sq_candidate;
1423  }
1424  }
1425  iterations++;
1426  }
1427 
1428  /* Set error if we didn't converge */
1429  if ( !( lambda < MAXLAMBDA && iterations < UVES_FIT_MAXITER ) )
1430  {
1431 #if DEBUG_LM
1432  uves_msg_error("Failed to converge, lambda=%f iterations=%d",
1433  lambda, iterations);
1434 #endif
1435  cpl_free(ia_local);
1436  cpl_matrix_delete(alpha);
1437  cpl_matrix_delete(beta);
1438  cpl_free(a_da);
1439  cpl_free(part);
1440  cpl_ensure_code(
1441  CPL_FALSE,
1442  CPL_ERROR_CONTINUE);
1443  }
1444 
1445  /* Compute mse if requested */
1446  if (mse != NULL)
1447  {
1448  int i;
1449 
1450  *mse = 0.0;
1451 
1452  for(i = 0; i < N; i++)
1453  {
1454  double fx_i = 0.0;
1455  double residual = 0.0;
1456 
1457  /* Evaluate f(x_i) at the best fit parameters */
1458  if( f(&(x_data[i*D]),
1459  a_data,
1460  &fx_i) != 0)
1461  {
1462  cpl_free(ia_local);
1463  cpl_matrix_delete(alpha);
1464  cpl_matrix_delete(beta);
1465  cpl_free(a_da);
1466  cpl_free(part);
1467  cpl_ensure_code(
1468  CPL_FALSE,
1469  CPL_ERROR_ILLEGAL_INPUT);
1470  }
1471 
1472  residual = y_data[i] - fx_i;
1473  *mse += residual * residual;
1474  }
1475  *mse /= N;
1476  }
1477 
1478  /* Compute reduced chi^2 if requested */
1479  if (red_chisq != NULL)
1480  {
1481  /* We already know the optimal chi^2 (and that N > Mfit)*/
1482  *red_chisq = chi_sq / (N-Mfit);
1483  }
1484 
1485  /* Compute covariance matrix if requested
1486  * cov = alpha(lambda=0)^-1
1487  */
1488  if (covariance != NULL)
1489  {
1490  cpl_matrix *cov;
1491 
1492  if( get_candidate(a_data, ia_local,
1493  M, N, D, 0.0, f, dfda,
1494  x_data, y_data, sigma_data,
1495  part, alpha, beta, a_da)
1496  != 0)
1497  {
1498  cpl_free(ia_local);
1499  cpl_matrix_delete(alpha);
1500  cpl_matrix_delete(beta);
1501  cpl_free(a_da);
1502  cpl_free(part);
1503  cpl_ensure_code(
1504  CPL_FALSE,
1505  cpl_error_get_code());
1506  }
1507 
1508  cov = cpl_matrix_invert_create(alpha);
1509  if (cov == NULL)
1510  {
1511  cpl_free(ia_local);
1512  cpl_matrix_delete(alpha);
1513  cpl_matrix_delete(beta);
1514  cpl_free(a_da);
1515  cpl_free(part);
1516  cpl_ensure_code(
1517  CPL_FALSE,
1518  cpl_error_get_code());
1519  }
1520 
1521  /* Make sure that variances are positive */
1522  {
1523  int i;
1524  for (i = 0; i < Mfit; i++)
1525  {
1526  if ( !(cpl_matrix_get(cov, i, i) > 0) )
1527  {
1528  cpl_free(ia_local);
1529  cpl_matrix_delete(alpha);
1530  cpl_matrix_delete(beta);
1531  cpl_free(a_da);
1532  cpl_free(part);
1533  cpl_matrix_delete(cov);
1534  *covariance = NULL;
1535  cpl_ensure_code(
1536  CPL_FALSE,
1537  CPL_ERROR_SINGULAR_MATRIX);
1538  }
1539  }
1540  }
1541 
1542  /* Expand covariance matrix from Mfit x Mfit
1543  to M x M. Set rows/columns corresponding to fixed
1544  parameters to zero */
1545 
1546  *covariance = cpl_matrix_new(M, M);
1547  if (*covariance == NULL)
1548  {
1549  cpl_free(ia_local);
1550  cpl_matrix_delete(alpha);
1551  cpl_matrix_delete(beta);
1552  cpl_free(a_da);
1553  cpl_free(part);
1554  cpl_matrix_delete(cov);
1555  cpl_ensure_code(
1556  CPL_FALSE,
1557  CPL_ERROR_ILLEGAL_OUTPUT);
1558  }
1559 
1560  {
1561  int j, jmfit;
1562 
1563  for (j = 0, jmfit = 0; j < M; j++)
1564  if (ia_local[j] != 0)
1565  {
1566  int i, imfit;
1567 
1568  for (i = 0, imfit = 0; i < M; i++)
1569  if (ia_local[i] != 0)
1570  {
1571  cpl_matrix_set(*covariance, i, j,
1572  cpl_matrix_get(
1573  cov, imfit, jmfit));
1574  imfit += 1;
1575  }
1576 
1577  assert( imfit == Mfit );
1578 
1579  jmfit += 1;
1580  }
1581 
1582  assert( jmfit == Mfit );
1583  }
1584 
1585  cpl_matrix_delete(cov);
1586  }
1587 
1588  cpl_free(ia_local);
1589  cpl_matrix_delete(alpha);
1590  cpl_matrix_delete(beta);
1591  cpl_free(a_da);
1592  cpl_free(part);
1593 
1594  return CPL_ERROR_NONE;
1595 }
1596 
1597 /*----------------------------------------------------------------------------*/
1607 /*----------------------------------------------------------------------------*/
1608 cpl_error_code
1609 uves_cast_image(cpl_image **image, cpl_type to_type)
1610 {
1611  cpl_image *temp = NULL;
1612 
1613  assure( image != NULL, CPL_ERROR_NULL_INPUT, "Null image");
1614 
1615  temp = *image;
1616 
1617  check( *image = cpl_image_cast(temp, to_type), "Couldn't convert image to %s",
1618  uves_tostring_cpl_type(to_type));
1619 
1620  cleanup:
1621  uves_free_image(&temp);
1622  return cpl_error_get_code();
1623 }
1624 
1625 
1626 /*-----------------------------------------------------------------------------*/
1641 /*-----------------------------------------------------------------------------*/
1642 cpl_error_code
1643 uves_crop_image(cpl_image **image, int x1, int y_1, int x2, int y2)
1644 {
1645  cpl_image *temp = NULL;
1646 
1647  assure( image != NULL, CPL_ERROR_NULL_INPUT, "Null image");
1648 
1649  temp = *image;
1650 
1651  check( *image = cpl_image_extract(temp, x1, y_1, x2, y2),
1652  "Could not extract image");
1653 
1654  cleanup:
1655  uves_free_image(&temp);
1656  return cpl_error_get_code();
1657 }
1658 
1659 /*----------------------------------------------------------------------------*/
1683 /*----------------------------------------------------------------------------*/
1684 cpl_error_code
1685 uves_get_property_value(const uves_propertylist *plist, const char *keyword,
1686  cpl_type keywordtype, void *result)
1687 {
1688  cpl_type t;
1689 
1690  /* Check input */
1691  assure( plist != NULL , CPL_ERROR_NULL_INPUT, "Null property list");
1692  assure( keyword != NULL, CPL_ERROR_NULL_INPUT, "Null keyword");
1693 
1694  /* Check for existence... */
1695  assure( uves_propertylist_contains(plist, keyword), CPL_ERROR_DATA_NOT_FOUND,
1696  "Keyword %s does not exist", keyword );
1697  /* ...and type of keyword */
1698  check( t = uves_propertylist_get_type(plist, keyword) ,
1699  "Could not read type of keyword '%s'", keyword );
1700  assure(t == keywordtype , CPL_ERROR_TYPE_MISMATCH,
1701  "Keyword '%s' has wrong type (%s). %s expected",
1702  keyword, uves_tostring_cpl_type(t), uves_tostring_cpl_type(keywordtype));
1703 
1704  /* Read the keyword */
1705  switch (keywordtype)
1706  {
1707  case CPL_TYPE_INT :
1708  check( *(( int *)result) =
1709  uves_propertylist_get_int(plist, keyword),
1710  "Could not get (integer) value of %s", keyword );
1711  break;
1712  case CPL_TYPE_BOOL :
1713  check( *(( bool *)result) =
1714  uves_propertylist_get_bool(plist, keyword),
1715  "Could not get (boolean) value of %s", keyword );
1716  break;
1717  case CPL_TYPE_DOUBLE:
1718  check( *(( double *)result) =
1719  uves_propertylist_get_double(plist, keyword),
1720  "Could not get (double) value of %s" , keyword );
1721  break;
1722  case CPL_TYPE_STRING:
1723  check( *((const char * *)result) =
1724  uves_propertylist_get_string(plist, keyword),
1725  "Could not get (string) value of %s" , keyword );
1726  break;
1727  default:
1728  assure( false, CPL_ERROR_INVALID_TYPE, "Unknown type");
1729  break;
1730  }
1731 
1732  cleanup:
1733  return cpl_error_get_code();
1734 }
1735 
1736 
1737 /*----------------------------------------------------------------------------*/
1762 /*----------------------------------------------------------------------------*/
1763 
1764 const char *
1765 uves_find_frame(const cpl_frameset *frames, const char **wanted, int N, int *found,
1766  const cpl_frame **frame)
1767 {
1768  const char *filename = NULL; /* Return NULL in case of error */
1769  int i;
1770  const cpl_frame *f = NULL;
1771 
1772  /* Return well-defined pointers in case of error */
1773  *found = 0;
1774  if (frame != NULL)
1775  {
1776  *frame = NULL;
1777  }
1778 
1779  for (i = 0; i < N; i++)
1780  {
1781  check( f = cpl_frameset_find_const(frames, wanted[i]),
1782  "Could not search through frame set");
1783  if (f != NULL)
1784  {
1785  check( filename = cpl_frame_get_filename(f),
1786  "Could not read frame filename");
1787  *found = i;
1788  if (frame != NULL)
1789  {
1790  *frame = f;
1791  }
1792  /* break */
1793  i = N;
1794  }
1795  }
1796 
1797  /* Set an error if a matching frame wasn't found */
1798  assure(filename != NULL, CPL_ERROR_DATA_NOT_FOUND, "No such frame in frame set");
1799 
1800  cleanup:
1801  return filename;
1802 }
1803 
1804 
1805 /*----------------------------------------------------------------------------*/
1812 /*----------------------------------------------------------------------------*/
1813 cpl_size
1814 uves_get_nextensions(const char *filename)
1815 {
1816  cpl_size result = 0;
1817  cpl_frame *f = NULL;
1818 
1819  /* CPL only supports reading the number of extensions of a FITS
1820  file if this is in a frame, so create a frame for the specified file */
1821 
1822  check(( f = cpl_frame_new(),
1823  cpl_frame_set_filename(f, filename)),
1824  "Could not create frame");
1825 
1826  check( result = cpl_frame_get_nextensions(f),
1827  "Error reading number of extensions of file '%s'", filename);
1828  cleanup:
1829  cpl_frame_delete(f);
1830  return result;
1831 }
1832 
1833 /*----------------------------------------------------------------------------*/
1854 /*----------------------------------------------------------------------------*/
1855 cpl_error_code
1856 uves_get_parameter(const cpl_parameterlist *parameters, const char *context,
1857  const char *recipe_id, const char *name, cpl_type type,
1858  void *value)
1859 {
1860  char *fullname = NULL;
1861  const cpl_parameter *p = NULL;
1862 
1863  passure( parameters != NULL, " ");
1864  /* 'context' may be NULL */
1865  passure( recipe_id != NULL, " ");
1866  passure( name != NULL, " ");
1867  passure( value != NULL, " ");
1868 
1869  if (context != NULL)
1870  {
1871  check( fullname = uves_sprintf("%s.%s.%s", context, recipe_id, name),
1872  "Error getting full parameter name");
1873  }
1874  else
1875  {
1876  check( fullname = uves_sprintf("%s.%s", recipe_id, name),
1877  "Error getting full parameter name");
1878  }
1879 
1880  /* Const cast */
1881  check( p = cpl_parameterlist_find_const(parameters, fullname),
1882  "Error searching for parameter '%s'", fullname);
1883  assure( p != NULL, CPL_ERROR_DATA_NOT_FOUND,
1884  "No parameter '%s' in parameter list", fullname);
1885 
1886  /* Check that parameter has the correct type */
1887  {
1888  cpl_type partype;
1889  check( partype = cpl_parameter_get_type(p),
1890  "Could not read type of parameter '%s'", fullname);
1891  assure( partype == type, CPL_ERROR_TYPE_MISMATCH,
1892  "Parameter '%s' has type %s. Expected type was %s",
1893  fullname,
1895  }
1896 
1897  /* Read the parameter */
1898  switch(type)
1899  {
1900  case CPL_TYPE_INT:
1901  check( *(int * )value =
1902  cpl_parameter_get_int (p),
1903  "Could not read integer parameter '%s'", fullname);
1904  break;
1905  case CPL_TYPE_BOOL:
1906  check( *(bool * )value =
1907  cpl_parameter_get_bool (p),
1908  "Could not read boolean parameter '%s'", fullname);
1909  break;
1910  case CPL_TYPE_DOUBLE:
1911  check( *(double * )value =
1912  cpl_parameter_get_double(p),
1913  "Could not read double parameter '%s'" , fullname );
1914  break;
1915  case CPL_TYPE_STRING:
1916  check( *(const char ** )value =
1917  cpl_parameter_get_string(p),
1918  "Could not read string parameter '%s'" , fullname );
1919  break;
1920  default:
1921  assure(false, CPL_ERROR_UNSUPPORTED_MODE,
1922  "Don't know how to read parameter '%s' of type %s",
1923  fullname, uves_tostring_cpl_type(type));
1924  break;
1925  }
1926 
1927  cleanup:
1928  cpl_free(fullname); fullname = NULL;
1929  return cpl_error_get_code();
1930 }
1931 
1932 /*----------------------------------------------------------------------------*/
1948 /*----------------------------------------------------------------------------*/
1949 cpl_error_code
1950 uves_set_parameter(cpl_parameterlist *parameters,
1951  const char *context,
1952  const char *name, cpl_type type, void *value)
1953 {
1954  char *fullname = NULL;
1955  cpl_parameter *p = NULL;
1956 
1957  check( fullname = uves_sprintf("%s.%s", context, name),
1958  "Error getting full parameter name");
1959 
1960  /* Const cast */
1961  check( p = cpl_parameterlist_find(parameters, fullname),
1962  "Error searching for parameter '%s'", fullname);
1963  assure( p != NULL, CPL_ERROR_DATA_NOT_FOUND,
1964  "No parameter '%s' in parameter list", fullname);
1965 
1966  /* Check that parameter has the correct type */
1967  {
1968  cpl_type partype;
1969  check( partype = cpl_parameter_get_type(p),
1970  "Could not read type of parameter '%s'", fullname);
1971  assure( partype == type, CPL_ERROR_TYPE_MISMATCH,
1972  "Parameter '%s' has type %s. Expected type was %s",
1973  fullname, uves_tostring_cpl_type(partype),
1974  uves_tostring_cpl_type(type));
1975  }
1976 
1977  /* Set the parameter */
1978  switch(type)
1979  {
1980  case CPL_TYPE_INT:
1981  check( cpl_parameter_set_int (p, *((int *) value)),
1982  "Could not set integer parameter '%s'", fullname);
1983  break;
1984  case CPL_TYPE_BOOL:
1985  check( cpl_parameter_set_bool (p, *((bool *) value)),
1986  "Could not set boolean parameter '%s'", fullname);
1987  break;
1988  case CPL_TYPE_DOUBLE:
1989  check( cpl_parameter_set_double(p, *((double *) value)),
1990  "Could not set double parameter '%s'" , fullname);
1991  break;
1992  case CPL_TYPE_STRING:
1993  check( cpl_parameter_set_string(p, *((char **) value)),
1994  "Could not set string parameter '%s'" , fullname);
1995  break;
1996  default:
1997  assure(false, CPL_ERROR_UNSUPPORTED_MODE,
1998  "Don't know how to set parameter of type %s",
1999  uves_tostring_cpl_type(type));
2000  break;
2001  }
2002 
2003  cleanup:
2004  cpl_free(fullname); fullname = NULL;
2005  return cpl_error_get_code();
2006 }
2007 
2008 /*----------------------------------------------------------------------------*/
2022 /*----------------------------------------------------------------------------*/
2023 
2024 cpl_error_code
2025 uves_set_parameter_default(cpl_parameterlist *parameters, const char *context,
2026  const char *parname,
2027  cpl_type type, void *value)
2028 {
2029  const char *full_name = NULL;
2030  cpl_parameter *p = NULL;
2031  cpl_type partype;
2032 
2033  if (context != NULL)
2034  {
2035  full_name = uves_sprintf("%s.%s", context, parname);
2036  }
2037  else
2038  {
2039  full_name = uves_sprintf("%s", parname);
2040  }
2041 
2042  if (full_name == NULL)
2043  {
2044  return CPL_ERROR_ILLEGAL_OUTPUT;
2045  }
2046 
2047 
2048  if ( (p = cpl_parameterlist_find(parameters, full_name)) == NULL)
2049  {
2050  uves_msg_error("Missing parameter: '%s'", full_name);
2051 
2052  uves_free_string_const(&full_name);
2053  if (cpl_error_get_code() != CPL_ERROR_NONE) return cpl_error_get_code();
2054  else return CPL_ERROR_DATA_NOT_FOUND;
2055  }
2056 
2057  partype = cpl_parameter_get_type(p);
2058 
2059  if (partype != type)
2060  {
2061  uves_msg_error("Parameter '%s' has type %s. Expected type was %s",
2062  full_name, uves_tostring_cpl_type(partype),
2063  uves_tostring_cpl_type(type));
2064 
2065  uves_free_string_const(&full_name);
2066  return CPL_ERROR_TYPE_MISMATCH;
2067  }
2068 
2069  switch(type)
2070  {
2071  case CPL_TYPE_INT:
2072  cpl_parameter_set_default_int (p, *((int *) value));
2073  break;
2074  case CPL_TYPE_BOOL:
2075  cpl_parameter_set_default_bool (p, *((bool *) value));
2076  break;
2077  case CPL_TYPE_DOUBLE:
2078  cpl_parameter_set_default_double(p, *((double *) value));
2079  break;
2080  case CPL_TYPE_STRING:
2081  cpl_parameter_set_default_string(p, *((char **) value));
2082  break;
2083  default:
2084  uves_msg_error("Unknown type: %s", uves_tostring_cpl_type(type));
2085 
2086  uves_free_string_const(&full_name);
2087  return CPL_ERROR_INVALID_TYPE;
2088  }
2089 
2090  if (cpl_error_get_code() != CPL_ERROR_NONE)
2091  {
2092  uves_msg_error("Error changing value of parameter '%s'",
2093  full_name);
2094 
2095  uves_free_string_const(&full_name);
2096  return cpl_error_get_code();
2097  }
2098 
2099 
2100  uves_free_string_const(&full_name);
2101  return CPL_ERROR_NONE;
2102 }
2103 
2104 
2105 /*----------------------------------------------------------------------------*/
2117 /*----------------------------------------------------------------------------*/
2118 void
2119 uves_raise_to_median_frac(cpl_table *t, const char *column, double fraction)
2120 {
2121  int i = 0;
2122  double threshold;
2123 
2124  assure_nomsg(t != NULL, CPL_ERROR_NULL_INPUT);
2125  assure(cpl_table_has_column(t, column), CPL_ERROR_DATA_NOT_FOUND,
2126  "No such column: %s", column);
2127 
2128  assure(cpl_table_get_column_type(t, column) == CPL_TYPE_DOUBLE,
2129  CPL_ERROR_UNSUPPORTED_MODE, "Column %s has type %s. %s expected",
2130  column,
2131  uves_tostring_cpl_type(cpl_table_get_column_type(t, column)),
2132  uves_tostring_cpl_type(CPL_TYPE_DOUBLE));
2133 
2134 
2135  threshold = fraction * cpl_table_get_column_median(t, column);
2136  for (i = 0; i < cpl_table_get_nrow(t); i++)
2137  {
2138  if (cpl_table_get_double(t, column, i, NULL) < threshold)
2139  {
2140  cpl_table_set_double(t, column, i, threshold);
2141  }
2142  }
2143 
2144  cleanup:
2145  return;
2146 }
2147 
2148 
2149 /*----------------------------------------------------------------------------*/
2166 /*----------------------------------------------------------------------------*/
2167 
2168 int
2169 uves_select_table_rows(cpl_table *t, const char *column,
2170  cpl_table_select_operator operator, double value)
2171 {
2172  int result = 0;
2173  cpl_type type;
2174 
2175  assure( t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
2176  assure( cpl_table_has_column(t, column), CPL_ERROR_INCOMPATIBLE_INPUT,
2177  "No such column: %s", column);
2178 
2179  type = cpl_table_get_column_type(t, column);
2180 
2181  assure( type == CPL_TYPE_DOUBLE || type == CPL_TYPE_FLOAT ||
2182  type == CPL_TYPE_INT, CPL_ERROR_INVALID_TYPE,
2183  "Column '%s' must be double or int. %s found", column,
2184  uves_tostring_cpl_type(type));
2185 
2186  check( cpl_table_select_all(t), "Error selecting rows");
2187 
2188  if (type == CPL_TYPE_DOUBLE)
2189  {
2190  result = cpl_table_and_selected_double(t, column, operator, value);
2191  }
2192  else if (type == CPL_TYPE_FLOAT)
2193  {
2194  result = cpl_table_and_selected_float(t, column, operator, value);
2195  }
2196  else if (type == CPL_TYPE_INT)
2197  {
2198  result = cpl_table_and_selected_int(t, column, operator,
2199  uves_round_double(value));
2200  }
2201  else { /*impossible*/ passure(false, ""); }
2202 
2203  cleanup:
2204  return result;
2205 
2206 }
2207 
2208 /*----------------------------------------------------------------------------*/
2222 /*----------------------------------------------------------------------------*/
2223 int
2224 uves_extract_table_rows_local(cpl_table *t, const char *column,
2225  cpl_table_select_operator operator, double value)
2226 {
2227  int result = 0;
2228 
2229  assure( t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
2230  assure( cpl_table_has_column(t, column), CPL_ERROR_INCOMPATIBLE_INPUT,
2231  "No such column: %s", column);
2232 
2233  check( result = uves_select_table_rows(t, column, operator, value),
2234  "Error selecting rows");
2235 
2236  check( cpl_table_not_selected(t), "Error selecting rows");
2237 
2238  check( cpl_table_erase_selected(t), "Error deleting rows");
2239 
2240  cleanup:
2241  return result;
2242 }
2243 
2244 /*----------------------------------------------------------------------------*/
2263 /*----------------------------------------------------------------------------*/
2264 int
2265 uves_erase_table_rows(cpl_table *t, const char *column,
2266  cpl_table_select_operator operator, double value)
2267 {
2268  int result = 0;
2269 
2270  assure( t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
2271  assure( cpl_table_has_column(t, column), CPL_ERROR_INCOMPATIBLE_INPUT,
2272  "No such column: %s", column);
2273 
2274  check( result = uves_select_table_rows(t, column, operator, value),
2275  "Error selecting rows");
2276 
2277  check( cpl_table_erase_selected(t), "Error deleting rows");
2278 
2279  cleanup:
2280  return result;
2281 }
2282 
2283 
2284 
2285 
2286 
2287 /*----------------------------------------------------------------------------*/
2296 /*----------------------------------------------------------------------------*/
2297 
2298 void uves_propertylist_append_property(uves_propertylist *plist, const cpl_property *p)
2299 {
2300  switch(cpl_property_get_type(p)) {
2301  case CPL_TYPE_CHAR:
2302  uves_propertylist_append_char(plist, cpl_property_get_name(p), cpl_property_get_char(p));
2303  break;
2304  case CPL_TYPE_BOOL:
2305  uves_propertylist_append_bool(plist, cpl_property_get_name(p), cpl_property_get_bool(p));
2306  break;
2307  case CPL_TYPE_INT:
2308  uves_propertylist_append_int(plist, cpl_property_get_name(p), cpl_property_get_int(p));
2309  break;
2310  case CPL_TYPE_LONG:
2311  uves_propertylist_append_long(plist, cpl_property_get_name(p), cpl_property_get_long(p));
2312  break;
2313  case CPL_TYPE_FLOAT:
2314  uves_propertylist_append_float(plist, cpl_property_get_name(p), cpl_property_get_float(p));
2315  break;
2316  case CPL_TYPE_DOUBLE:
2317  uves_propertylist_append_double(plist, cpl_property_get_name(p), cpl_property_get_double(p));
2318  break;
2319  case CPL_TYPE_STRING:
2320  uves_propertylist_append_string(plist, cpl_property_get_name(p), cpl_property_get_string(p));
2321  break;
2322  default:
2323  assure( false, CPL_ERROR_UNSUPPORTED_MODE,
2324  "Type is %s", uves_tostring_cpl_type(cpl_property_get_type(p)));
2325  break;
2326  }
2327  cleanup:
2328  return;
2329 }
2330 
2331 
2332 
2333 /*----------------------------------------------------------------------------*/
2339 /*----------------------------------------------------------------------------*/
2340 int
2341 uves_table_and_selected_invalid(cpl_table *t, const char *column)
2342 {
2343  if (cpl_table_get_column_type(t, column) != CPL_TYPE_STRING)
2344  {
2345  return cpl_table_and_selected_invalid(t, column);
2346  }
2347  else
2348  {
2349  int i = 0;
2350  for (i = 0; i < cpl_table_get_nrow(t); i++)
2351  {
2352  if (cpl_table_is_selected(t, i))
2353  {
2354  if (cpl_table_is_valid(t, column, i))
2355  {
2356  cpl_table_unselect_row(t, i);
2357  }
2358  /* else keep it selected */
2359  }
2360  /* else unselected, don't change */
2361  }
2362 
2363  return cpl_table_count_selected(t);
2364  }
2365 }
2366 
2367 
2368 /*----------------------------------------------------------------------------*/
2385 /*----------------------------------------------------------------------------*/
2386 int
2387 uves_erase_invalid_table_rows(cpl_table *t, const char *column)
2388 {
2389  int result = 0;
2390 
2391  assure( t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
2392 
2393  if (column == NULL)
2394  /* Loop through all columns */
2395  {
2396  const char *name;
2397  cpl_table *argument = t;
2398  result = 0;
2399  while ( (name = cpl_table_get_column_name(argument)) != NULL)
2400  {
2401  int n_deleted_rows;
2402 
2403  argument = NULL;
2404  n_deleted_rows = uves_erase_invalid_table_rows(t, name);
2405 
2406  if (n_deleted_rows > 0)
2407  {
2408  uves_msg_low("%d rows with invalid '%s' removed",
2409  n_deleted_rows, name);
2410  }
2411  result += n_deleted_rows;
2412  }
2413  }
2414  else
2415  {
2416  assure( cpl_table_has_column(t, column), CPL_ERROR_INCOMPATIBLE_INPUT,
2417  "No such column: %s", column);
2418 
2419  check(( cpl_table_select_all(t),
2420  result = uves_table_and_selected_invalid(t, column), /* workaround here */
2421  cpl_table_erase_selected(t)), /* and here */
2422  "Error deleting rows");
2423  }
2424 
2425  cleanup:
2426  return result;
2427 }
2428 
2429 /*----------------------------------------------------------------------------*/
2446 /*----------------------------------------------------------------------------*/
2447 cpl_table *
2448 uves_extract_table_rows(const cpl_table *t, const char *column,
2449  cpl_table_select_operator operator, double value)
2450 {
2451  cpl_table *result = NULL;
2452 
2453  assure( t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
2454  assure( cpl_table_has_column(t, column), CPL_ERROR_INCOMPATIBLE_INPUT,
2455  "No such column: %s", column);
2456 
2457  /* 1. Extract (duplicate) the entire table
2458  2. remove rows *not* satisfying the criterion */
2459  check(( result = cpl_table_duplicate(t),
2460 
2461  uves_select_table_rows(result, column, operator, value),
2462  cpl_table_not_selected(result), /* Inverses selection */
2463 
2464  cpl_table_erase_selected(result)),
2465 
2466  "Error extracting rows");
2467 
2468  passure( cpl_table_count_selected(result) == cpl_table_get_nrow(result),
2469  "%" CPL_SIZE_FORMAT " %" CPL_SIZE_FORMAT "",
2470  cpl_table_count_selected(result), cpl_table_get_nrow(result) );
2471 
2472  cleanup:
2473  if (cpl_error_get_code() != CPL_ERROR_NONE)
2474  {
2475  uves_free_table(&result);
2476  }
2477  return result;
2478 }
2479 
2480 /*----------------------------------------------------------------------------*/
2492 /*----------------------------------------------------------------------------*/
2493 void
2494 uves_sort_table_1(cpl_table *t, const char *column, bool reverse)
2495 {
2496  uves_propertylist *plist = NULL;
2497 
2498  assure(t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
2499  assure(cpl_table_has_column(t, column), CPL_ERROR_ILLEGAL_INPUT,
2500  "No column '%s'", column);
2501 
2502  check(( plist = uves_propertylist_new(),
2503  uves_propertylist_append_bool(plist, column, reverse)),
2504  "Could not create property list for sorting");
2505 
2506  check( uves_table_sort(t, plist), "Could not sort table");
2507 
2508  cleanup:
2509  uves_free_propertylist(&plist);
2510  return;
2511 }
2512 
2513 /*----------------------------------------------------------------------------*/
2529 /*----------------------------------------------------------------------------*/
2530 void
2531 uves_sort_table_2(cpl_table *t, const char *column1, const char *column2,
2532  bool reverse1, bool reverse2)
2533 {
2534  uves_propertylist *plist = NULL;
2535 
2536  assure(t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
2537  assure(cpl_table_has_column(t, column1), CPL_ERROR_ILLEGAL_INPUT,
2538  "No column '%s'", column1);
2539  assure(cpl_table_has_column(t, column2), CPL_ERROR_ILLEGAL_INPUT,
2540  "No column '%s'", column2);
2541 
2542  check(( plist = uves_propertylist_new(),
2543  uves_propertylist_append_bool(plist, column1, reverse1),
2544  uves_propertylist_append_bool(plist, column2, reverse2)),
2545  "Could not create property list for sorting");
2546  check( uves_table_sort(t, plist), "Could not sort table");
2547 
2548  cleanup:
2549  uves_free_propertylist(&plist);
2550  return;
2551 }
2552 /*----------------------------------------------------------------------------*/
2567 /*----------------------------------------------------------------------------*/
2568 void
2569 uves_sort_table_3(cpl_table *t, const char *column1, const char *column2,
2570  const char *column3,
2571  bool reverse1, bool reverse2, bool reverse3)
2572 {
2573  uves_propertylist *plist = NULL;
2574 
2575  assure(t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
2576  assure(cpl_table_has_column(t, column1), CPL_ERROR_ILLEGAL_INPUT,
2577  "No column '%s'", column1);
2578  assure(cpl_table_has_column(t, column2), CPL_ERROR_ILLEGAL_INPUT,
2579  "No column '%s'", column2);
2580  assure(cpl_table_has_column(t, column3), CPL_ERROR_ILLEGAL_INPUT,
2581  "No column '%s'", column3);
2582 
2583  check(( plist = uves_propertylist_new(),
2584  uves_propertylist_append_bool(plist, column1, reverse1),
2585  uves_propertylist_append_bool(plist, column2, reverse2),
2586  uves_propertylist_append_bool(plist, column3, reverse3)),
2587  "Could not create property list for sorting");
2588  check( uves_table_sort(t, plist), "Could not sort table");
2589 
2590  cleanup:
2591  uves_free_propertylist(&plist);
2592  return;
2593 }
2594 /*----------------------------------------------------------------------------*/
2599 /*----------------------------------------------------------------------------*/
2600 void uves_free(const void *mem)
2601 {
2602  cpl_free((void *)mem); /* No, it is not a bug. The cast is safe */
2603  return;
2604 }
2605 
2606 /*----------------------------------------------------------------------------*/
2611 /*----------------------------------------------------------------------------*/
2612 void uves_free_image(cpl_image **i)
2613 {if(i){cpl_image_delete(*i); *i = NULL;}}
2614 /*----------------------------------------------------------------------------*/
2619 /*----------------------------------------------------------------------------*/
2620 void uves_free_mask(cpl_mask **m)
2621 {if(m){cpl_mask_delete(*m); *m = NULL;}}
2622 /*----------------------------------------------------------------------------*/
2627 /*----------------------------------------------------------------------------*/
2628 void uves_free_imagelist(cpl_imagelist **i)
2629 {if(i){cpl_imagelist_delete(*i); *i = NULL;}}
2630 
2631 /*----------------------------------------------------------------------------*/
2636 /*----------------------------------------------------------------------------*/
2637 void uves_free_table(cpl_table **t)
2638 {if(t){cpl_table_delete(*t); *t = NULL;}}
2639 
2640 /*----------------------------------------------------------------------------*/
2645 /*----------------------------------------------------------------------------*/
2646 void uves_free_table_const(const cpl_table **t)
2647 {if(t){cpl_table_delete((cpl_table*) (*t)); *t = NULL;}}
2648 
2649 /*----------------------------------------------------------------------------*/
2654 /*----------------------------------------------------------------------------*/
2655 void uves_free_propertylist(uves_propertylist **p)
2656 {if(p){uves_propertylist_delete(*p); *p = NULL;}}
2657 
2658 /*----------------------------------------------------------------------------*/
2663 /*----------------------------------------------------------------------------*/
2664 void uves_free_propertylist_const(const uves_propertylist **p)
2665 {if(p){uves_propertylist_delete(*p); *p = NULL;}}
2666 
2667 /*----------------------------------------------------------------------------*/
2672 /*----------------------------------------------------------------------------*/
2673 void uves_free_property(cpl_property **p)
2674 {if(p){cpl_property_delete(*p); *p = NULL;}}
2675 /*----------------------------------------------------------------------------*/
2680 /*----------------------------------------------------------------------------*/
2681 void uves_free_polynomial(cpl_polynomial **p)
2682 {if(p){cpl_polynomial_delete(*p); *p = NULL;}}
2683 /*----------------------------------------------------------------------------*/
2688 /*----------------------------------------------------------------------------*/
2689 void uves_free_matrix(cpl_matrix **m)
2690 {if(m){cpl_matrix_delete(*m); *m = NULL;}}
2691 /*----------------------------------------------------------------------------*/
2696 /*----------------------------------------------------------------------------*/
2697 void uves_free_parameterlist(cpl_parameterlist **p)
2698 {if(p){cpl_parameterlist_delete(*p); *p = NULL;}}
2699 /*----------------------------------------------------------------------------*/
2704 /*----------------------------------------------------------------------------*/
2705 void uves_free_frameset(cpl_frameset **f)
2706 {if(f){cpl_frameset_delete(*f); *f = NULL;}}
2707 /*----------------------------------------------------------------------------*/
2712 /*----------------------------------------------------------------------------*/
2713 void uves_free_frame(cpl_frame **f)
2714 {if(f){cpl_frame_delete(*f); *f = NULL;}}
2715 /*----------------------------------------------------------------------------*/
2720 /*----------------------------------------------------------------------------*/
2721 void uves_free_bivector(cpl_bivector **b)
2722 {if(b){cpl_bivector_delete(*b); *b = NULL;}}
2723 /*----------------------------------------------------------------------------*/
2728 /*----------------------------------------------------------------------------*/
2729 void uves_free_vector(cpl_vector **v)
2730 {if(v){cpl_vector_delete(*v); *v = NULL;}}
2731 /*----------------------------------------------------------------------------*/
2736 /*----------------------------------------------------------------------------*/
2737 void uves_free_stats(cpl_stats **s)
2738 {if(s){cpl_stats_delete(*s); *s = NULL;}}
2739 /*----------------------------------------------------------------------------*/
2744 /*----------------------------------------------------------------------------*/
2745 void uves_unwrap_vector(cpl_vector **v)
2746 {if(v){cpl_vector_unwrap(*v); *v = NULL;}}
2747 /*----------------------------------------------------------------------------*/
2752 /*----------------------------------------------------------------------------*/
2753 void uves_unwrap_vector_const(const cpl_vector **v)
2754 {if(v){cpl_vector_unwrap((cpl_vector*) (*v)); *v = NULL;}}
2755 /*----------------------------------------------------------------------------*/
2760 /*----------------------------------------------------------------------------*/
2761 void uves_unwrap_bivector_vectors(cpl_bivector **b)
2762 {if(b){cpl_bivector_unwrap_vectors(*b); *b = NULL;}}
2763 /*----------------------------------------------------------------------------*/
2768 /*----------------------------------------------------------------------------*/
2769 void uves_free_array(cpl_array **a)
2770 {if(a){cpl_array_delete(*a); *a = NULL;}}
2771 
2772 /*----------------------------------------------------------------------------*/
2777 /*----------------------------------------------------------------------------*/
2778 void uves_free_int(int **i)
2779 {if(i){cpl_free(*i); *i = NULL;}}
2780 
2781 /*----------------------------------------------------------------------------*/
2786 /*----------------------------------------------------------------------------*/
2787 void uves_free_int_const(const int **i)
2788 {if(i){uves_free(*i); *i = NULL;}}
2789 /*----------------------------------------------------------------------------*/
2794 /*----------------------------------------------------------------------------*/
2795 void uves_free_float(float **f)
2796 {if(f){cpl_free(*f); *f = NULL;}}
2797 
2798 /*----------------------------------------------------------------------------*/
2803 /*----------------------------------------------------------------------------*/
2804 void uves_free_double(double **d)
2805 {if(d){cpl_free(*d); *d = NULL;}}
2806 
2807 /*----------------------------------------------------------------------------*/
2812 /*----------------------------------------------------------------------------*/
2813 void uves_free_string(char **s)
2814 {if(s){cpl_free(*s); *s = NULL;}}
2815 /*----------------------------------------------------------------------------*/
2820 /*----------------------------------------------------------------------------*/
2821 void uves_free_string_const(const char **s)
2822 {if(s){uves_free(*s); *s = NULL;}}
2823 
2824 /*----------------------------------------------------------------------------*/
2845 /*----------------------------------------------------------------------------*/
2846 
2847 static double
2848 get_chisq(int N, int D,
2849  int (*f)(const double x[], const double a[], double *result),
2850  const double *a,
2851  const double *x,
2852  const double *y,
2853  const double *sigma)
2854 {
2855  double chi_sq; /* Result */
2856  int i = 0;
2857 
2858  /* For efficiency, don't check input in this static function */
2859 
2860  chi_sq = 0.0;
2861  for (i = 0; i < N; i++)
2862  {
2863  double fx_i;
2864  double residual; /* Residual in units of uncertainty */
2865  const double *x_i = &(x[0+i*D]);
2866 
2867  /* Evaluate */
2868  cpl_ensure( f(x_i,
2869  a,
2870  &fx_i) == 0, CPL_ERROR_ILLEGAL_INPUT, -1.0);
2871 
2872  /* Accumulate */
2873  if (sigma == NULL)
2874  {
2875  residual = (fx_i - y[i]);
2876  }
2877  else
2878  {
2879  residual = (fx_i - y[i]) / sigma[i];
2880  }
2881 
2882  chi_sq += residual*residual;
2883 
2884  }
2885 
2886  return chi_sq;
2887 }
2888 /*----------------------------------------------------------------------------*/
2922 /*----------------------------------------------------------------------------*/
2923 static int
2924 get_candidate(const double *a, const int ia[],
2925  int M, int N, int D,
2926  double lambda,
2927  int (*f)(const double x[], const double a[], double *result),
2928  int (*dfda)(const double x[], const double a[], double result[]),
2929  const double *x,
2930  const double *y,
2931  const double *sigma,
2932  double *partials,
2933  cpl_matrix *alpha,
2934  cpl_matrix *beta,
2935  double *a_da)
2936 {
2937  int Mfit = 0; /* Number of non-constant fit parameters */
2938  cpl_matrix *da; /* Solution of alpha * da = beta */
2939  double *alpha_data;
2940  double *beta_data;
2941  double *da_data;
2942  int i, imfit = 0;
2943  int j, jmfit = 0;
2944  int k = 0;
2945 
2946  /* For efficiency, don't check input in this static function */
2947 
2948  Mfit = cpl_matrix_get_nrow(alpha);
2949 
2950  alpha_data = cpl_matrix_get_data(alpha);
2951  beta_data = cpl_matrix_get_data(beta);
2952 
2953  /* Build alpha, beta:
2954  *
2955  * alpha[i,j] = sum_{k=1,N} (sigma_k)^-2 * df/da_i * df/da_j *
2956  * (1 + delta_ij lambda) ,
2957  *
2958  * beta[i] = sum_{k=1,N} (sigma_k)^-2 * ( y_k - f(x_k) ) * df/da_i
2959  *
2960  * where (i,j) loop over the non-constant parameters (0 to Mfit-1),
2961  * delta is Kronecker's delta, and all df/da are evaluated in x_k
2962  */
2963 
2964  cpl_matrix_fill(alpha, 0.0);
2965  cpl_matrix_fill(beta , 0.0);
2966 
2967  for (k = 0; k < N; k++)
2968  {
2969  double sm2 = 0.0; /* (sigma_k)^-2 */
2970  double fx_k = 0.0; /* f(x_k) */
2971  const double *x_k = &(x[0+k*D]); /* x_k */
2972 
2973  if (sigma == NULL)
2974  {
2975  sm2 = 1.0;
2976  }
2977  else
2978  {
2979  sm2 = 1.0 / (sigma[k] * sigma[k]);
2980  }
2981 
2982  /* Evaluate f(x_k) */
2983  cpl_ensure( f(x_k, a, &fx_k) == 0, CPL_ERROR_ILLEGAL_INPUT, -1);
2984 
2985  /* Evaluate (all) df/da (x_k) */
2986  cpl_ensure( dfda(x_k, a, partials) == 0,
2987  CPL_ERROR_ILLEGAL_INPUT, -1);
2988 
2989  for (i = 0, imfit = 0; i < M; i++)
2990  {
2991  if (ia[i] != 0)
2992  {
2993  /* Beta */
2994  beta_data[imfit] +=
2995  sm2 * (y[k] - fx_k) * partials[i];
2996 
2997  /* Alpha is symmetrical, so compute
2998  only lower-left part */
2999  for (j = 0, jmfit = 0; j < i; j++)
3000  {
3001  if (ia[j] != 0)
3002  {
3003  alpha_data[jmfit + imfit*Mfit] +=
3004  sm2 * partials[i] *
3005  partials[j];
3006 
3007  jmfit += 1;
3008  }
3009  }
3010 
3011  /* Alpha, diagonal terms */
3012  j = i;
3013  jmfit = imfit;
3014 
3015  alpha_data[jmfit + imfit*Mfit] +=
3016  sm2 * partials[i] *
3017  partials[j] * (1 + lambda);
3018 
3019  imfit += 1;
3020  }
3021  }
3022 
3023  assert( imfit == Mfit );
3024  }
3025 
3026  /* Create upper-right part of alpha */
3027  for (i = 0, imfit = 0; i < M; i++)
3028  {
3029  if (ia[i] != 0)
3030  {
3031  for (j = i+1, jmfit = imfit+1; j < M; j++)
3032  {
3033  if (ia[j] != 0)
3034  {
3035  alpha_data[jmfit + imfit*Mfit] =
3036  alpha_data[imfit + jmfit*Mfit];
3037 
3038  jmfit += 1;
3039  }
3040  }
3041  assert( jmfit == Mfit );
3042 
3043  imfit += 1;
3044  }
3045  }
3046  assert( imfit == Mfit );
3047 
3048  da = cpl_matrix_solve(alpha, beta);
3049 
3050  cpl_ensure(da != NULL, cpl_error_get_code(), -1);
3051 
3052  /* Create a+da vector by adding a and da */
3053  da_data = cpl_matrix_get_data(da);
3054 
3055  for (i = 0, imfit = 0; i < M; i++)
3056  {
3057  if (ia[i] != 0)
3058  {
3059  a_da[i] = a[i] + da_data[0 + imfit*1];
3060 
3061  imfit += 1;
3062  }
3063  else
3064  {
3065  a_da[i] = a[i];
3066  }
3067  }
3068 
3069  assert( imfit == Mfit );
3070 
3071  cpl_matrix_delete(da);
3072 
3073  return 0;
3074 }
3075