SINFONI Pipeline Reference Manual  2.6.0
sinfo_recipes.c
1 /*
2  * This file is part of the ESO SINFONI Pipeline
3  * Copyright (C) 2004,2005 European Southern Observatory
4  *
5  * This program 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 * E.S.O. - VLT project
21 *
22 *
23 *
24 * who when what
25 * -------- -------- ----------------------------------------------
26 * schreib 05/06/00 created
27 */
28 
29 #ifdef HAVE_CONFIG_H
30 # include <config.h>
31 #endif
32 
33 #include "sinfo_vltPort.h"
34 
35 /*
36  * System Headers
37  */
38 
39 /*
40  * Local Headers
41  */
42 
43 #include "sinfo_recipes.h"
44 #include "sinfo_globals.h"
45 
46 /*----------------------------------------------------------------------------
47  * Local variables
48  *--------------------------------------------------------------------------*/
49 
50 static float sqrarg ;
51 
52 static double chi1 ; /* old reduced chi-squared */
53 static double sinfo_chi2 ; /* new reduced chi-squared */
54 static double labda ; /* mixing parameter */
55 static double vec[MAXPAR] ; /* correction sinfo_vector */
56 static double matrix1[MAXPAR][MAXPAR] ; /* original sinfo_matrix */
57 static double matrix2[MAXPAR][MAXPAR] ; /* inverse of matrix1 */
58 static int nfree ; /* number of free parameters */
59 static int parptr[MAXPAR] ; /* parameter pointer */
60 
61 /*----------------------------------------------------------------------------
62  * Defines
63  *--------------------------------------------------------------------------*/
64 
65 #define SQR(a) (sqrarg = (a) , sqrarg*sqrarg)
66 
67 /*----------------------------------------------------------------------------
68  * Functions private to this module
69  *--------------------------------------------------------------------------*/
70 
71 
72 static int new_inv_mat (void) ;
73 
74 static void new_get_mat ( float * xdat,
75  int * xdim,
76  float * ydat,
77  float * wdat,
78  int * ndat,
79  float * fpar,
80  float * epar/*,
81  int * npar */) ;
82 
83 static int new_get_vec ( float * xdat,
84  int * xdim,
85  float * ydat,
86  float * wdat,
87  int * ndat,
88  float * fpar,
89  float * epar,
90  int * npar ) ;
91 
92 static float
93 new_gaussian ( float * xdat, float * parlist/*, int * npar*/ );
94 static void
95 new_gaussian_deriv( float * xdat, float * parlist,
96  float * dervs/*, int * npar*/ );
97 
98 
99 
107 /*----------------------------------------------------------------------------
108  * Function codes
109  *--------------------------------------------------------------------------*/
110 
111 
112 float sinfo_new_f_median(float * array, int n)
113 {
114  pixelvalue p_array[100];
115  int i;
116 
117  for (i=0;i<n;i++)
118  p_array[i]= (pixelvalue) array[i];
119 
120  return (float) sinfo_new_median(p_array, n);
121 }
122 
123 
140 float sinfo_new_clean_mean( float * array,
141  int n_elements,
142  float throwaway_low,
143  float throwaway_high )
144 {
145  int i, n ;
146  int lo_n, hi_n ;
147  float sum ;
148 
149  if ( array == NULL )
150  {
151  sinfo_msg_error(" no array given in sinfo_clean_mean!") ;
152  return FLT_MAX ;
153  }
154 
155  if ( n_elements <= 0 )
156  {
157  sinfo_msg_error("wrong number of elements given") ;
158  return FLT_MAX ;
159  }
160 
161  if ( throwaway_low < 0. || throwaway_high < 0. ||
162  throwaway_low + throwaway_high >= 100. )
163  {
164  sinfo_msg_error("wrong throw away percentage given!") ;
165  return FLT_MAX ;
166  }
167 
168  lo_n = (int) (throwaway_low * (float)n_elements / 100.) ;
169  hi_n = (int) (throwaway_high * (float)n_elements / 100.) ;
170 
171  /* sort the array */
172  sinfo_pixel_qsort( array, n_elements ) ;
173 
174  n = 0 ;
175  sum = 0. ;
176  for ( i = lo_n ; i < n_elements - hi_n ; i++ )
177  {
178  if ( !isnan(array[i]) )
179  {
180  sum += array[i] ;
181  n++ ;
182  }
183  }
184  if ( n == 0 )
185  {
186  return FLAG ;
187  }
188  else
189  {
190  return sum/(float)n ;
191  }
192 }
193 
194 /*--------------------------------------------------------------------------*/
207 pixelvalue sinfo_new_median(pixelvalue * array, int n)
208 {
209  pixelvalue med ;
210 
211  if ( array == NULL || n <= 0 )
212  {
213  sinfo_msg_warning("nothing in the pixelvalue array, ZERO returned");
214  return ZERO ;
215  }
216 
217  if ( n == 1 )
218  {
219  return array[0] ;
220  }
221 
222  sinfo_pixel_qsort((float*) array, n) ;
223  if ( n % 2 == 1 )
224  {
225  med = array[n/2] ;
226  }
227  else
228  {
229  med = (array[n/2] + array[n/2 - 1])/2. ;
230  }
231  return med ;
232 }
233 
234 
235 
236 
237 
285 int sinfo_new_lsqfit_c ( float * xdat,
286  int * xdim,
287  float * ydat,
288  float * wdat,
289  int * ndat,
290  float * fpar,
291  float * epar,
292  int * mpar,
293  int * npar,
294  float * tol ,
295  int * its ,
296  float * lab )
297 {
298  int i, n, r ;
299  int itc ; /* fate of fit */
300  int found ; /* fit converged: 1, not yet converged: 0 */
301  int nuse ; /* number of useable data points */
302  double tolerance ; /* accuracy */
303 
304  itc = 0 ; /* fate of fit */
305  found = 0 ; /* reset */
306  nfree = 0 ; /* number of free parameters */
307  nuse = 0 ; /* number of legal data points */
308 
309  if ( *tol < (FLT_EPSILON * 10.0 ) )
310  {
311  tolerance = FLT_EPSILON * 10.0 ; /* default tolerance */
312  }
313  else
314  {
315  tolerance = *tol ; /* tolerance */
316  }
317 
318  labda = fabs( *lab ) * LABFAC ; /* start value for mixing parameter */
319  for ( i = 0 ; i < (*npar) ; i++ )
320  {
321  if ( mpar[i] )
322  {
323  if ( nfree > MAXPAR ) /* too many free parameters */
324  {
325  return -1 ;
326  }
327  parptr[nfree++] = i ; /* a free parameter */
328  }
329  }
330 
331  if (nfree == 0) /* no free parameters */
332  {
333  return -2 ;
334  }
335 
336  for ( n = 0 ; n < (*ndat) ; n++ )
337  {
338  if ( wdat[n] > 0.0 ) /* legal weight */
339  {
340  nuse ++ ;
341  }
342  }
343 
344  if ( nfree >= nuse )
345  {
346  return -3 ; /* no degrees of freedom */
347  }
348  if ( labda == 0.0 ) /* linear fit */
349  {
350  /* initialize fpar array */
351  for ( i = 0 ; i < nfree ; fpar[parptr[i++]] = 0.0 ) ;
352  new_get_mat ( xdat, xdim, ydat, wdat, ndat, fpar, epar/*, npar*/ ) ;
353  r = new_get_vec ( xdat, xdim, ydat, wdat, ndat, fpar, epar, npar ) ;
354  if ( r ) /* error */
355  {
356  return r ;
357  }
358  for ( i = 0 ; i < (*npar) ; i++ )
359  {
360  fpar[i] = epar[i] ; /* save new parameters */
361  epar[i] = 0.0 ; /* and set errors to zero */
362  }
363  chi1 = sqrt( chi1 / (double) (nuse - nfree) ) ;
364  for ( i = 0 ; i < nfree ; i++ )
365  {
366  if ( (matrix1[i][i] <= 0.0 ) || (matrix2[i][i] <= 0.0) )
367  {
368  return -7 ;
369  }
370  epar[parptr[i]] = chi1 * sqrt( matrix2[i][i] ) /
371  sqrt( matrix1[i][i] ) ;
372  }
373  }
374  else /* non-linear fit */
375  {
376  /*----------------------------------------------------------------
377  * the non-linear fit uses the steepest descent method in combination
378  * with the Taylor method. The mixing of these methods is controlled
379  * by labda. In the outer loop ( called the iteration loop ) we build
380  * the matrix and calculate the correction vector. In the inner loop
381  * (called the interpolation loop) we check whether we have obtained a
382  * better solution than the previous one. If so, we leave the inner
383  loop else we increase labda (give more weight to the steepest
384  descent method) calculate the correction vector and check again.
385  After the inner loop we do a final check on the goodness of the
386  fit and if this satisfies
387  * the tolerance we calculate the errors of the fitted parameters.
388  */
389  while ( !found ) /* iteration loop */
390  {
391  if ( itc++ == (*its) ) /* increase iteration counter */
392  {
393  return -4 ;
394  }
395  new_get_mat( xdat, xdim, ydat, wdat, ndat, fpar, epar/*, npar */) ;
396 
397  /*-------------------------------------------------------------
398  * here we decrease labda since we may assume that each iteration
399  * brings us closer to the answer.
400  */
401  if ( labda > LABMIN )
402  {
403  labda = labda / LABFAC ; /* decrease labda */
404  }
405  r = new_get_vec( xdat, xdim, ydat, wdat, ndat, fpar, epar, npar ) ;
406  if ( r ) /* error */
407  {
408  return r ;
409  }
410 
411  while ( chi1 >= sinfo_chi2 ) /* interpolation loop */
412  {
413  /*-----------------------------------------------------------
414  * The next statement is based on experience, not on the
415  mathematics of the problem. It is assumed that we have
416  reached convergence when the pure steepest descent method
417  does not produce a better solution.
418  */
419  if ( labda > LABMAX ) /* assume solution found */
420  {
421  break ;
422  }
423  labda = labda * LABFAC ; /* increase mixing parameter */
424  r = new_get_vec(xdat,xdim,ydat,wdat,ndat,fpar,epar,npar) ;
425  if ( r ) /* error */
426  {
427  return r ;
428  }
429  }
430 
431  if ( labda <= LABMAX ) /* save old parameters */
432  {
433  for ( i = 0 ; i < *npar ; i++ )
434  {
435  fpar[i] = epar[i] ;
436  }
437  }
438  if ( (fabs( sinfo_chi2 - chi1 ) <= (tolerance * chi1)) ||
439  (labda > LABMAX) )
440  {
441  /*-------------------------------------------------------------
442  * we have a satisfying solution, so now we need to calculate
443  the correct errors of the fitted parameters. This we do by
444  using the pure Taylor method because we are very close to
445  the real solution.
446  */
447  labda = 0.0 ; /* for Taylor solution */
448  new_get_mat(xdat,xdim,ydat,wdat,ndat,fpar,epar/*, npar */) ;
449  r = new_get_vec(xdat,xdim,ydat,wdat,ndat,fpar,epar,npar) ;
450 
451  if ( r ) /* error */
452  {
453  return r ;
454  }
455  for ( i = 0 ; i < (*npar) ; i++ )
456  {
457  epar[i] = 0.0 ; /* set error to zero */
458  }
459  sinfo_chi2 = sqrt ( sinfo_chi2 / (double) (nuse - nfree) ) ;
460 
461  for ( i = 0 ; i < nfree ; i++ )
462  {
463  if ( (matrix1[i][i] <= 0.0) || (matrix2[i][i] <= 0.0) )
464  {
465  return -7 ;
466  }
467  epar[parptr[i]] = sinfo_chi2 * sqrt( matrix2[i][i] ) /
468  sqrt( matrix1[i][i] ) ;
469  }
470  found = 1 ; /* we found a solution */
471  }
472  }
473  }
474  return itc ; /* return number of iterations */
475 }
476 
477 
478 
485 void sinfo_new_convert_ZEROs_to_0_for_images(cpl_image * im)
486 {
487  int i ;
488  int ilx=0;
489  int ily=0;
490  float* pidata=NULL;
491 
492  if ( im == NULL )
493  {
494  sinfo_msg_error("no input image given!\n") ;
495  return ;
496  }
497  ilx=cpl_image_get_size_x(im);
498  ily=cpl_image_get_size_y(im);
499  pidata=cpl_image_get_data(im);
500  for ( i = 0 ; i < (int) ilx*ily ; i++ )
501  {
502  if( isnan(pidata[i]) )
503  {
504  pidata[i] = 0. ;
505  }
506  }
507  return ;
508 }
509 
518 void sinfo_new_convert_ZEROs_to_0_for_cubes(cpl_imagelist * cube)
519 {
520  int i ;
521  int inp=0;
522 
523  if ( cube == NULL )
524  {
525  sinfo_msg_error("no input cube given!") ;
526  return ;
527  }
528  inp=cpl_imagelist_get_size(cube);
529 
530  for ( i = 0 ; i < inp ; i++ )
531  {
532  cpl_image* i_img=cpl_imagelist_get(cube,i);
533  sinfo_new_convert_ZEROs_to_0_for_images(i_img) ;
534  cpl_imagelist_set(cube,i_img,i);
535  }
536  return ;
537 }
538 
539 
548 void
549 sinfo_new_convert_ZEROs_to_0_for_cubes_range(cpl_imagelist * cube,
550  const int z_min,const int z_max)
551 {
552  int i ;
553 
554  if ( cube == NULL )
555  {
556  sinfo_msg_error("no input cube given!") ;
557  return ;
558  }
559  for ( i = z_min ; i < z_max ; i++ )
560  {
561  cpl_image* i_img=cpl_imagelist_get(cube,i);
562  sinfo_new_convert_ZEROs_to_0_for_images(i_img) ;
563  cpl_imagelist_set(cube,i_img,i);
564  }
565  return ;
566 }
573 void sinfo_new_convert_0_to_ZEROs_for_images(cpl_image * im)
574 {
575  int i ;
576  int ilx=0;
577  int ily=0;
578  float* pidata=NULL;
579 
580  if ( im == NULL )
581  {
582  sinfo_msg_error("no input image given!") ;
583  return ;
584  }
585  ilx=cpl_image_get_size_x(im);
586  ily=cpl_image_get_size_y(im);
587  pidata=cpl_image_get_data(im);
588  for ( i = 0 ; i < (int) ilx*ily ; i++ )
589  {
590  if( pidata[i] == 0. )
591  {
592  pidata[i] = ZERO ;
593  }
594  }
595  return ;
596 }
597 
604 void sinfo_new_convert_0_to_ZERO_for_cubes(cpl_imagelist * cube)
605 {
606  int i ;
607  int inp=0;
608 
609  if ( cube == NULL )
610  {
611  sinfo_msg_error("no input cube given!") ;
612  return ;
613  }
614  inp=cpl_imagelist_get_size(cube);
615  for ( i = 0 ; i < inp ; i++ )
616  {
617  cpl_image* i_img=cpl_imagelist_get(cube,i);
618  sinfo_new_convert_0_to_ZEROs_for_images(i_img) ;
619  cpl_imagelist_set(cube,i_img,i);
620  }
621  return ;
622 }
623 
624 
633 void
634 sinfo_new_convert_0_to_ZERO_for_cubes_range(cpl_imagelist * cube,
635  const int z_min,const int z_max)
636 {
637  int i ;
638  //int inp=0;
639 
640  if ( cube == NULL )
641  {
642  sinfo_msg_error("no input cube given!") ;
643  return ;
644  }
645  //inp=cpl_imagelist_get_size(cube);
646  for ( i = z_min ; i < z_max ; i++ )
647  {
648  cpl_image* i_img=cpl_imagelist_get(cube,i);
649  sinfo_new_convert_0_to_ZEROs_for_images(i_img) ;
650  cpl_imagelist_set(cube,i_img,i);
651  }
652  return ;
653 }
660 void sinfo_new_invert(cpl_image * im)
661 {
662  int i ;
663  int ilx=0;
664  int ily=0;
665  float* pidata=NULL;
666 
667  ilx=cpl_image_get_size_x(im);
668  ily=cpl_image_get_size_y(im);
669  pidata=cpl_image_get_data(im);
670 
671  for ( i = 0 ; i < (int) ilx*ily ; i++ )
672  {
673  pidata[i] = -pidata[i] ;
674  }
675  return ;
676 }
677 
685 int sinfo_new_nint ( double x )
686 {
687  int k ;
688 
689  k = x ;
690  if ( x >= 0. )
691  {
692  if ( (x - (double) k) <= 0.5 )
693  {
694  return k ;
695  }
696  else
697  {
698  return k + 1 ;
699  }
700  }
701  else
702  {
703  if ( (x - (double) k) <= -0.5 )
704  {
705  return k - 1;
706  }
707  else
708  {
709  return k ;
710  }
711  }
712 }
713 
714 
728 #define STEP_MIN (-half_search)
729 #define STEP_MAX (half_search)
730 
731 double * sinfo_new_xcorrel(
732  float * line_i,
733  int width_i,
734  float * line_t,
735  int width_t,
736  int half_search,
737  int * delta,
738  int * maxpos,
739  double * xcorr_max
740 
741 )
742 {
743  double * xcorr ;
744  double mean_i, mean_t ;
745  double rms_i, rms_t ;
746  double sum, sqsum ;
747  double norm ;
748  int nsteps ;
749  int i ;
750  int step ;
751  /*double r;*/
752 
753 
754  /* Compute normalization factors */
755  sum = sqsum = 0.00 ;
756  for (i=0 ; i<width_i ; i++) {
757  sum += (double)line_i[i] ;
758  sqsum += (double)line_i[i] * (double)line_i[i];
759  }
760  mean_i = sum / (double)width_i ;
761  sqsum /= (double)width_i ;
762  rms_i = sqsum - mean_i*mean_i ;
763 
764  sum = sqsum = 0.00 ;
765  for (i=0 ; i<width_t ; i++) {
766  sum += (double)line_t[i] ;
767  sqsum += (double)line_t[i] * (double)line_t[i];
768  }
769  mean_t = sum / (double)width_t ;
770  sqsum /= (double)width_t ;
771  rms_t = sqsum - mean_t*mean_t ;
772 
773  norm = 1.00 / sqrt(rms_i * rms_t);
774 
775  nsteps = (STEP_MAX - STEP_MIN) ;
776  xcorr = cpl_malloc(nsteps * sizeof(double));
777  for (step=STEP_MIN ; step<STEP_MAX ; step++) {
778  xcorr[step-STEP_MIN] = 0.00 ;
779  int nval = 0 ;
780  for (i=0 ; i<width_t ; i++) {
781  if ((i+step >= 0) &&
782  (i+step < width_i)) {
783  xcorr[step-STEP_MIN] += ((double)line_t[i] - mean_t) *
784  ((double)line_i[i+step] - mean_i) *
785  norm ;
786  nval++ ;
787  }
788  }
789  xcorr[step-STEP_MIN] /= (double)nval ;
790  }
791  *xcorr_max = xcorr[0] ;
792  *maxpos = 0 ;
793  for (i=0 ; i<nsteps ; i++) {
794  if (xcorr[i]>*xcorr_max) {
795  *maxpos = i ;
796  *xcorr_max = xcorr[i];
797  }
798  }
799  (*delta) = + (STEP_MIN + *maxpos);
800  return xcorr ;
801 }
802 
803 /* FILE ELEMENT: sinfo_nev_ille.c */
804 /* */
805 /**********************************************************************/
806 /* */
807 /* double sinfo_nev_ille() */
808 /* */
809 /**********************************************************************/
810 /* */
811 /* DESCRIPTION: */
812 /* For a given table (x , f(x )), i = 0(1)n and a given argument z */
813 /* the function computes the interpolated value for the argument z */
814 /* using Neville's interpolation/ extrapolation algorithm. */
815 /* */
816 /* FUNCTIONS CALLED: */
817 /* System library: <stdio.h> printf(), fabs(); */
818 /* Numlib library: None */
819 /* Local functions: nevtable(); */
820 /* User supplied: None */
821 /* */
822 /* PROGRAMMED BY: T.Haavie */
823 /* DATE/VERSION: 88-07-06/1.0 */
824 /* */
825 /**********************************************************************/
826 double sinfo_nev_ille(double x[], double f[], int n, double z, int* flag)
827  /* PARAMETERS(input): */
828 /* double x[]; Abscissae values in the table. */
829 /* double f[]; Function values in the table. */
830 /* int n; The number of elements in the table is n+1. */
831 /* double z; Argument to be used in interpolation/extrapolation. */
832 
833 
834 /* PARAMETERS(input/output): */
835 /* int *flag; Flag parameter(output): */
836  /* = 0, n < 0 and/or eps < 0, should be positive+. */
837  /* = 1, required rel.err. is not obtained. */
838  /* = 2, required rel. err. is obtained. */
839 
840 /* the computed estimate for the interpolated/extrapolated value re- */
841 /* turned through function name sinfo_nev_ille. */
842 
843 {
844  double p[11]; /* Array used for storing the new row elements */
845  /* in the interpolation/extrapolation table. */
846  double q[11]; /* Array used for storing the old row elements */
847  /* in the interpolation/extrapolation table */
848 
849  double factor;
850  int m, k;
851 
852 
853 
854  if (n < 0 )
855  {
856  *flag = 0;
857  return(0.);
858  }
859 
860 
861  q[0] = f[0]; /* Set initial value in the table. */
862 
863  for (k = 1; k <= n; k++) /* k counts rows in the table. */
864  {
865  p[0] = f[k];
866  for (m = 1; m <= k; m++) /* m counts element in row. */
867  {
868  factor = (z - x[k]) / (x[k] - x[k-m]);
869  p[m] = p[m-1] + factor * (p[m-1] - q[m-1]);
870  }
871 
872 
873  for (m = 0; m <= k; m++) /* Shift old row to new row. */
874  q[m] = p[m];
875 
876  } /* End of k-loop. */
877 
878  *flag = 1; /* Required rel.error is not obtained. */
879  return(p[n]);
880 
881 } /* End of sinfo_nev_ille(). */
882 
883 
884 
885 float sinfo_new_nev_ille(float x[], float f[], int n, float z, int* flag)
886  /* PARAMETERS(input): */
887 /* float x[]; Abscissae values in the table. */
888 /* float f[]; Function values in the table. */
889 /* int n; The number of elements in the table is n+1. */
890 /* float z; Argument to be used in interpolation/extrapolation. */
891 
892 
893 /* PARAMETERS(input/output): */
894 /* int *flag; Flag parameter(output): */
895  /* = 0, n < 0 and/or eps < 0, should be positive+. */
896  /* = 1, required rel.err. is not obtained. */
897  /* = 2, required rel. err. is obtained. */
898 
899 /* the computed estimate for the interpolated/extrapolated value re- */
900 /* turned through function name sinfo_nev_ille. */
901 
902 {
903  float p[11]; /* Array used for storing the new row elements */
904  /* in the interpolation/extrapolation table. */
905  float q[11]; /* Array used for storing the old row elements */
906  /* in the interpolation/extrapolation table */
907 
908  float factor;
909  int m, k;
910 
911 
912 
913  if (n < 0 )
914  {
915  *flag = 0;
916  return(0.);
917  }
918 
919 
920  q[0] = f[0]; /* Set initial value in the table. */
921 
922  for (k = 1; k <= n; k++) /* k counts rows in the table. */
923  {
924  p[0] = f[k];
925  for (m = 1; m <= k; m++) /* m counts element in row. */
926  {
927  factor = (z - x[k]) / (x[k] - x[k-m]);
928  p[m] = p[m-1] + factor * (p[m-1] - q[m-1]);
929  }
930 
931 
932  for (m = 0; m <= k; m++) /* Shift old row to new row. */
933  q[m] = p[m];
934 
935  } /* End of k-loop. */
936 
937  *flag = 1; /* Required rel.error is not obtained. */
938  return(p[n]);
939 
940 } /* End of sinfo_nev_ille(). */
941 
942 
967 static int new_get_vec ( float * xdat,
968  int * xdim,
969  float * ydat,
970  float * wdat,
971  int * ndat,
972  float * fpar,
973  float * epar,
974  int * npar )
975 {
976  double dy ;
977  double mii ;
978  double mji ;
979  double mjj ;
980  int i, j, n, r ;
981 
982  /* loop to modify and scale the sinfo_matrix */
983  for ( j = 0 ; j < nfree ; j++ )
984  {
985  mjj = matrix1[j][j] ;
986  if ( mjj <= 0.0 ) /* diagonal element wrong */
987  {
988  return -5 ;
989  }
990  mjj = sqrt( mjj ) ;
991  for ( i = 0 ; i < j ; i++ )
992  {
993  mji = matrix1[j][i] / mjj / sqrt( matrix1[i][i] ) ;
994  matrix2[i][j] = matrix2[j][i] = mji ;
995  }
996  matrix2[j][j] = 1.0 + labda ; /* scaled value on diagonal */
997  }
998 
999  if ( (r = new_inv_mat()) ) /* sinfo_invert sinfo_matrix inlace */
1000  {
1001  return r ;
1002  }
1003 
1004  for ( i = 0 ; i < (*npar) ; i ++ )
1005  {
1006  epar[i] = fpar[i] ;
1007  }
1008 
1009  /* loop to calculate correction sinfo_vector */
1010  for ( j = 0 ; j < nfree ; j++ )
1011  {
1012  double dj = 0.0 ;
1013  mjj = matrix1[j][j] ;
1014  if ( mjj <= 0.0) /* not allowed */
1015  {
1016  return -7 ;
1017  }
1018  mjj = sqrt ( mjj ) ;
1019  for ( i = 0 ; i < nfree ; i++ )
1020  {
1021  mii = matrix1[i][i] ;
1022  if ( mii <= 0.0 )
1023  {
1024  return -7 ;
1025  }
1026  mii = sqrt( mii ) ;
1027  dj += vec[i] * matrix2[j][i] / mjj / mii ;
1028  }
1029  epar[parptr[j]] += dj ; /* new parameters */
1030  }
1031  chi1 = 0.0 ; /* reset reduced chi-squared */
1032 
1033  /* loop through the data points */
1034  for ( n = 0 ; n < (*ndat) ; n++ )
1035  {
1036  double wn = wdat[n] ; /* get weight */
1037  if ( wn > 0.0 ) /* legal weight */
1038  {
1039  dy = ydat[n] - new_gaussian( &xdat[(*xdim) * n], epar/*, npar*/ ) ;
1040  chi1 += wdat[n] * dy * dy ;
1041  }
1042  }
1043  return 0 ;
1044 }
1045 
1046 
1062 static void new_get_mat ( float * xdat,
1063  int * xdim,
1064  float * ydat,
1065  float * wdat,
1066  int * ndat,
1067  float * fpar,
1068  float * epar/*,
1069  int * npar */)
1070 {
1071 
1072  for ( int j = 0 ; j < nfree ; j++ )
1073  {
1074  vec[j] = 0.0 ; /* zero sinfo_vector */
1075  for ( int i = 0 ; i<= j ; i++ )
1076  /* zero sinfo_matrix only on and below diagonal */
1077  {
1078  matrix1[j][i] = 0.0 ;
1079  }
1080  }
1081  sinfo_chi2 = 0.0 ; /* reset reduced chi-squared */
1082 
1083  /* loop through data points */
1084  for ( int n = 0 ; n < (*ndat) ; n++ )
1085  {
1086  double wn = wdat[n] ;
1087  if ( wn > 0.0 ) /* legal weight ? */
1088  {
1089  double yd = ydat[n] - new_gaussian( &xdat[(*xdim) * n], fpar/*, npar*/ ) ;
1090  new_gaussian_deriv( &xdat[(*xdim) * n], fpar, epar/*, npar*/ ) ;
1091  sinfo_chi2 += yd * yd * wn ; /* add to chi-squared */
1092  for ( int j = 0 ; j < nfree ; j++ )
1093  {
1094  double wd = epar[parptr[j]] * wn ; /* weighted derivative */
1095  vec[j] += yd * wd ; /* fill sinfo_vector */
1096  for ( int i = 0 ; i <= j ; i++ ) /* fill sinfo_matrix */
1097  {
1098  matrix1[j][i] += epar[parptr[i]] * wd ;
1099  }
1100  }
1101  }
1102  }
1103 }
1104 
1105 
1106 
1107 
1108 
1109 
1118 static int new_inv_mat (void)
1119 {
1120  double even ;
1121  double hv[MAXPAR] ;
1122  double mjk ;
1123  int evin ;
1124  int i, j, k ;
1125  int per[MAXPAR] ;
1126 
1127  /* set permutation array */
1128  for ( i = 0 ; i < nfree ; i++ )
1129  {
1130  per[i] = i ;
1131  }
1132 
1133  for ( j = 0 ; j < nfree ; j++ ) /* in j-th column */
1134  {
1135  /* determine largest element of a row */
1136  double rowmax = fabs ( matrix2[j][j] ) ;
1137  int row = j ;
1138 
1139  for ( i = j + 1 ; i < nfree ; i++ )
1140  {
1141  if ( fabs ( matrix2[i][j] ) > rowmax )
1142  {
1143  rowmax = fabs( matrix2[i][j] ) ;
1144  row = i ;
1145  }
1146  }
1147 
1148  /* determinant is zero! */
1149  if ( matrix2[row][j] == 0.0 )
1150  {
1151  return -6 ;
1152  }
1153 
1154  /*if the largest element is not on the diagonal, then permutate rows */
1155  if ( row > j )
1156  {
1157  for ( k = 0 ; k < nfree ; k++ )
1158  {
1159  even = matrix2[j][k] ;
1160  matrix2[j][k] = matrix2[row][k] ;
1161  matrix2[row][k] = even ;
1162  }
1163  /* keep track of permutation */
1164  evin = per[j] ;
1165  per[j] = per[row] ;
1166  per[row] = evin ;
1167  }
1168 
1169  /* modify column */
1170  even = 1.0 / matrix2[j][j] ;
1171  for ( i = 0 ; i < nfree ; i++ )
1172  {
1173  matrix2[i][j] *= even ;
1174  }
1175  matrix2[j][j] = even ;
1176 
1177  for ( k = 0 ; k < j ; k++ )
1178  {
1179  mjk = matrix2[j][k] ;
1180  for ( i = 0 ; i < j ; i++ )
1181  {
1182  matrix2[i][k] -= matrix2[i][j] * mjk ;
1183  }
1184  for ( i = j + 1 ; i < nfree ; i++ )
1185  {
1186  matrix2[i][k] -= matrix2[i][j] * mjk ;
1187  }
1188  matrix2[j][k] = -even * mjk ;
1189  }
1190 
1191  for ( k = j + 1 ; k < nfree ; k++ )
1192  {
1193  mjk = matrix2[j][k] ;
1194  for ( i = 0 ; i < j ; i++ )
1195  {
1196  matrix2[i][k] -= matrix2[i][j] * mjk ;
1197  }
1198  for ( i = j + 1 ; i < nfree ; i++ )
1199  {
1200  matrix2[i][k] -= matrix2[i][j] * mjk ;
1201  }
1202  matrix2[j][k] = -even * mjk ;
1203  }
1204  }
1205 
1206  /* finally, repermute the columns */
1207  for ( i = 0 ; i < nfree ; i++ )
1208  {
1209  for ( k = 0 ; k < nfree ; k++ )
1210  {
1211  hv[per[k]] = matrix2[i][k] ;
1212  }
1213  for ( k = 0 ; k < nfree ; k++ )
1214  {
1215  matrix2[i][k] = hv[k] ;
1216  }
1217  }
1218 
1219  /* all is well */
1220  return 0 ;
1221 }
1222 
1223 
1224 
1225 
1246 float new_gaussian ( float * xdat, float * parlist/*, int * npar*/ )
1247 {
1248  double xd ; /* FWHM's of gauss function */
1249  double x ; /* position */
1250 
1251  xd = fabs((double) parlist[1]) ;
1252  x = (double) xdat[0] - (double) parlist[2] ;
1253  return (float) (
1254  (double) parlist[0] * exp( -4.0 * log(2.0) * (x/xd) * (x/xd) )
1255  + (double) parlist[3] ) ;
1256 }
1257 
1258 
1283 void
1284 new_gaussian_deriv(float * xdat,float * parlist,float * dervs/*, int * npar*/ )
1285 {
1286  double xd ; /* FWHM of sinfo_gaussian */
1287  double x, expon ; /* position and exponent */
1288 
1289  xd = fabs( (double) parlist[1] ) ;
1290 
1291  /* offset from peak position */
1292  x = (double) xdat[0] - (double) parlist[2] ;
1293 
1294  /* determine the derivatives: */
1295  expon = -4.0 * log(2.0) * (x/xd) * (x/xd) ;
1296  expon = exp( expon ) ;
1297 
1298  /* partial derivative by the amplitude */
1299  dervs[0] = (float) expon ;
1300 
1301  /* calculate a * exp(-arg) */
1302  expon = (double) parlist[0] * expon ;
1303 
1304  /* partial derivative by FWHM */
1305  dervs[1] = (float) ( expon * 8.0 * log(2.0) * x*x / (xd*xd*xd) ) ;
1306 
1307  /* partial derivative by the x position (parlist[2]) */
1308  dervs[2] = (float) (expon * 8.0 * log(2.0) * x/(xd*xd) ) ;
1309 
1310  /* partial derivative by the zero level */
1311  dervs[3] = 1.0 ;
1312 }
1313 
1314 
1315 /*==================================================================*/
1316 
1317 
1337 void
1338 sinfo_my_fit(float x[], float y[], int ndata, float sig[], int mwt, float *a,
1339  float *b, float *siga, float *sigb, float *chi2, float *q)
1340 {
1341  int i ;
1342  float t, sxoss, sx=0., sy=0., st2=0., ss ;
1343 
1344  *b = 0. ; /*accumulate sums ...*/
1345  if ( mwt )
1346  {
1347  ss = 0. ;
1348  for ( i = 0 ; i < ndata ; i++ ) /*... with weights*/
1349  {
1350  float wt = 1./SQR(sig[i]) ;
1351  ss += wt ;
1352  sx += x[i]*wt ;
1353  sy += y[i]*wt ;
1354  }
1355  }
1356  else
1357  {
1358  for ( i = 0 ; i < ndata ; i++ ) /*... or without weights*/
1359  {
1360  sx += x[i] ;
1361  sy += y[i] ;
1362  }
1363  ss = ndata ;
1364  }
1365  sxoss = sx/ss ;
1366 
1367  if ( mwt )
1368  {
1369  for ( i = 0 ; i < ndata ; i ++ )
1370  {
1371  t = (x[i] - sxoss)/sig[i] ;
1372  st2 += t*t ;
1373  *b += t*y[i]/sig[i] ;
1374  }
1375  }
1376  else
1377  {
1378  for ( i = 0 ; i < ndata ; i++ )
1379  {
1380  t = x[i] - sxoss ;
1381  st2 += t*t ;
1382  *b += t*y[i] ;
1383  }
1384  }
1385 
1386  *b /= st2 ;
1387  *a = (sy - sx*(*b))/ss ;
1388  *siga = sqrt ((1.0 + sx*sx/(ss*st2))/ss) ;
1389  *sigb = sqrt (1.0/st2) ;
1390  *chi2 = 0.0 ; /*calculate chi-square*/
1391  if ( mwt == 0 )
1392  {
1393  for ( i = 0 ; i < ndata ; i++ )
1394  {
1395  *chi2 += SQR (y[i] - (*a) - (*b)*x[i]) ;
1396  }
1397  *q = 1. ;
1398 
1399  /*------------------------------------------------------------------
1400  * for unweighted data evaluate typical sig using chi2, and adjust
1401  * the standard deviation
1402  */
1403  float sigdat = sqrt ((*chi2)/(ndata - 2)) ;
1404  *siga *= sigdat ;
1405  *sigb *= sigdat ;
1406  }
1407  else
1408  {
1409  for (i = 0 ; i < ndata ; i++)
1410  {
1411  *chi2 += SQR ((y[i] - (*a) - (*b) * x[i])/sig[i]) ;
1412  }
1413  *q = 1. ; /* delete rest of lines. q is not a good value */
1414  }
1415 }
1416 
1431 int sinfo_new_correlation ( float * data1, float * data2, int ndata )
1432 {
1433  /*float help[3*ndata] ;
1434  float corsum[3*ndata] ;*/
1435  float* help=NULL ;
1436  float* corsum=NULL ;
1437  float maxval ;
1438  int i, j, k, position, shift ;
1439  int /*start,end,size,*/ndata3,limit;
1440 
1441 
1442  /*ndata3=3*ndata;*/
1443  ndata3=ndata+300;
1444 
1445  if ( NULL == data1 || NULL == data2 || ndata <= 1 )
1446  {
1447  sinfo_msg_error(" wrong input for sinfo_correlation\n") ;
1448  return INT32_MAX ;
1449  }
1450 
1451  /* initialize the help arrays with zeros */
1452  help=cpl_calloc(ndata+300,sizeof(float));
1453  for ( i = 0 ; i < ndata3 ; i++ )
1454  {
1455  help[i] = 0. ;
1456  }
1457 
1458  /* shift the second data array by ndata in the help array */
1459  for ( i = 0 ; i < ndata ; i++ )
1460  {
1461  help[(300/2) + i] = data2[i] ;
1462  }
1463 
1464  /* compute the cross sinfo_correlation sum array */
1465  corsum=cpl_calloc(ndata+300,sizeof(float));
1466  for ( j = 0 ; j < ndata3 ; j++ )
1467  {
1468  if ( ndata3-j <= ndata)
1469  limit = ndata3-j;
1470  else
1471  limit = ndata;
1472  corsum[j] = 0. ;
1473  for ( k = 0 ; k < limit ; k++ )
1474  {
1475  /*if ( k + j >= ndata3 )
1476  {
1477  break ;
1478  }*/
1479  corsum[j] += data1[k] * help[k + j] ;
1480  }
1481  }
1482 
1483  /* search for the maximal corsum value and determine its position */
1484  maxval = -FLT_MAX ;
1485  position = -1 ;
1486  for ( i = 0 ; i < ndata3 ; i++ )
1487  {
1488  if ( maxval < corsum[i] )
1489  {
1490  maxval = corsum[i] ;
1491  position = i ;
1492  }
1493  }
1494 
1495  /* determine shift of data2 relative to the data1 array */
1496  shift = position - 300/2 ;
1497  cpl_free(help);
1498  cpl_free(corsum);
1499 
1500  return shift ;
1501 }
1502 
1503 /*--------------------------------------------------------------------------*/
#define sinfo_msg_error(...)
Print an error message.
Definition: sinfo_msg.h:69
#define sinfo_msg_warning(...)
Print an warning message.
Definition: sinfo_msg.h:93