SINFONI Pipeline Reference Manual  2.6.0
sinfo_focus.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 16/01/02 created
27 */
28 
29 /************************************************************************
30 * NAME
31 * sinfo_focus.c -
32 *------------------------------------------------------------------------
33 */
34 
35 #ifdef HAVE_CONFIG_H
36 # include <config.h>
37 #endif
38 #include "sinfo_vltPort.h"
39 
40 /*
41  * System Headers
42  */
43 
44 /*
45  * Local Headers
46  */
47 
48 #include "sinfo_focus.h"
49 #include "sinfo_recipes.h"
50 #include <float.h>
51 
52 /*----------------------------------------------------------------------------
53  * Defines
54  *--------------------------------------------------------------------------*/
55 
56 #define XDIMG 2 /* dimension of the x values */
57 #define TOLG 0.001 /* fitting tolerance */
58 #define LABG 0.1 /* labda parameter */
59 #define ITSG 200 /* maximum number of iterations */
60 #define LABFACG 10.0 /* labda step factor */
61 #define LABMAXG 1.0e+10 /* maximum value for labda */
62 #define LABMING 1.0e-10 /* minimum value for labda */
63 #define NPAR 7 /* number of fit parameters */
64 #define PI_NUMB (3.1415926535897932384626433832795) /* pi */
65 
66 
67 /*----------------------------------------------------------------------------
68  * Local variables
69  *--------------------------------------------------------------------------*/
70 
71 static double chi1 ; /* old reduced chi-squared */
72 static double chi2 ; /* new reduced chi-squared */
73 static double labda ; /* mixing parameter */
74 static double vec[NPAR] ; /* correction sinfo_vector */
75 static double matrix1[NPAR][NPAR] ; /* original sinfo_matrix */
76 static double matrix2[NPAR][NPAR] ; /* inverse of matrix1 */
77 static int nfree ; /* number of free parameters */
78 static int parptr[NPAR] ; /* parameter pointer */
79 
80 /*----------------------------------------------------------------------------
81  * Functions private to this module
82  *--------------------------------------------------------------------------*/
83 
84 static int new_inv_mat (void) ;
85 
86 static void new_get_mat ( double * xdat,
87  int * xdim,
88  double * ydat,
89  double * wdat,
90  int * ndat,
91  double * fpar,
92  double * epar/*,
93  int * npar */) ;
94 
95 static int new_get_vec ( double * xdat,
96  int * xdim,
97  double * ydat,
98  double * wdat,
99  int * ndat,
100  double * fpar,
101  double * epar,
102  int * npar ) ;
103 
104 static int new_gauss2Ellipse ( double * parlist ) ;
113 /*----------------------------------------------------------------------------
114  * Function codes
115  *--------------------------------------------------------------------------*/
116 
117 /*-------------------------------------------------------------------------*/
141 /*--------------------------------------------------------------------------*/
142 
143 double sinfo_new_gaussian_ellipse(double * xdat, double * parlist)
144 {
145  double result ;
146  double x ;
147  double y ;
148  double fwhmx ;
149  double fwhmy ;
150  double costheta ;
151  double sintheta ;
152  double argX ; /* arguments in the exponent */
153  double argY ;
154 
155  /* some abbreviations */
156  x = xdat[0] - parlist[0] ;
157  y = xdat[1] - parlist[1] ;
158 
159  fwhmx = fabs(parlist[4]) ;
160  fwhmy = fabs(parlist[5]) ;
161 
162  costheta = cos ( parlist[6] ) ;
163  sintheta = sin ( parlist[6] ) ;
164 
165  argX = x * costheta + y * sintheta ;
166  argY = -x * sintheta + y * costheta ;
167 
168  /* function */
169  result = parlist[2] * exp(-4.*log(2.0)*((argX/fwhmx)*(argX/fwhmx)+
170  (argY/fwhmy)*(argY/fwhmy))) +
171  parlist[3] ;
172 
173  return result ;
174 }
175 
200 void
201 sinfo_new_gaussian_ellipse_deriv(double * xdat,
202  double * parlist,
203  double * dervs )
204 {
205  double x ;
206  double y ;
207  double fwhmx ;
208  double fwhmy ;
209  double argX ;
210  double argY ;
211  double expon ;
212  double e8log2 ;
213  double fwx2 ;
214  double fwy2 ;
215  double costheta ;
216  double sintheta ;
217 
218  /* some abbreviations */
219  x = xdat[0] - parlist[0] ;
220  y = xdat[1] - parlist[1] ;
221 
222  fwhmx = fabs(parlist[4]) ;
223  fwhmy = fabs(parlist[5]) ;
224  fwx2 = fwhmx * fwhmx ;
225  fwy2 = fwhmy * fwhmy ;
226 
227  costheta = cos ( parlist[6] ) ;
228  sintheta = sin ( parlist[6] ) ;
229 
230  argX = x * costheta + y * sintheta ;
231  argY = -x * sintheta + y * costheta ;
232 
233  expon = exp ( -4.0 * log(2.0) * ((argX/fwhmx)*(argX/fwhmx) +
234  (argY/fwhmy)*(argY/fwhmy)) ) ;
235  e8log2 = expon * 8.0 * log(2.0) ;
236 
237  /* determine the derivatives */
238  /* partial derivative x-position */
239  dervs[0] = -parlist[2]*e8log2 * (-argX*costheta/fwx2 + argY*sintheta/fwy2);
240  /* partial derivative y-position */
241  dervs[1] = -parlist[2]*e8log2 * (-argX*sintheta/fwx2 - argY*costheta/fwy2);
242  /* partial derivative amplitude */
243  dervs[2] = expon ;
244  /* partial derivative background */
245  dervs[3] = 1. ;
246  /* partial derivative fwhmx */
247  dervs[4] = parlist[2]*e8log2 * argX*argX/(fwx2*fwhmx) ;
248  /* partial derivative fwhmy */
249  dervs[5] = parlist[2]*e8log2 * argY*argY/(fwy2*fwhmy) ;
250  /* partial derivative theta */
251  dervs[6] = -parlist[2]*e8log2 * argY * argX * (1.0/fwx2 - 1.0/fwy2) ;
252 
253 }
254 
266 static int new_inv_mat (void)
267 {
268  double even ;
269  double hv[NPAR] ;
270  double mjk ;
271 
272  int evin ;
273  int i, j, k;
274  int per[NPAR] ;
275 
276  /* set permutation array */
277  for ( i = 0 ; i < nfree ; i++ )
278  {
279  per[i] = i ;
280  }
281 
282  for ( j = 0 ; j < nfree ; j++ ) /* in j-th column */
283  {
284  /* determine largest element of a row */
285  double rowmax = fabs ( matrix2[j][j] ) ;
286  int row = j ;
287 
288  for ( i = j + 1 ; i < nfree ; i++ )
289  {
290  if ( fabs ( matrix2[i][j] ) > rowmax )
291  {
292  rowmax = fabs( matrix2[i][j] ) ;
293  row = i ;
294  }
295  }
296 
297  /* determinant is zero! */
298  if ( matrix2[row][j] == 0.0 )
299  {
300  return -6 ;
301  }
302 
303  /* if the largest element is not on the diagonal,
304  then permutate rows */
305  if ( row > j )
306  {
307  for ( k = 0 ; k < nfree ; k++ )
308  {
309  even = matrix2[j][k] ;
310  matrix2[j][k] = matrix2[row][k] ;
311  matrix2[row][k] = even ;
312  }
313  /* keep track of permutation */
314  evin = per[j] ;
315  per[j] = per[row] ;
316  per[row] = evin ;
317  }
318 
319  /* modify column */
320  even = 1.0 / matrix2[j][j] ;
321  for ( i = 0 ; i < nfree ; i++ )
322  {
323  matrix2[i][j] *= even ;
324  }
325  matrix2[j][j] = even ;
326 
327  for ( k = 0 ; k < j ; k++ )
328  {
329  mjk = matrix2[j][k] ;
330  for ( i = 0 ; i < j ; i++ )
331  {
332  matrix2[i][k] -= matrix2[i][j] * mjk ;
333  }
334  for ( i = j + 1 ; i < nfree ; i++ )
335  {
336  matrix2[i][k] -= matrix2[i][j] * mjk ;
337  }
338  matrix2[j][k] = -even * mjk ;
339  }
340 
341  for ( k = j + 1 ; k < nfree ; k++ )
342  {
343  mjk = matrix2[j][k] ;
344  for ( i = 0 ; i < j ; i++ )
345  {
346  matrix2[i][k] -= matrix2[i][j] * mjk ;
347  }
348  for ( i = j + 1 ; i < nfree ; i++ )
349  {
350  matrix2[i][k] -= matrix2[i][j] * mjk ;
351  }
352  matrix2[j][k] = -even * mjk ;
353  }
354  }
355 
356  /* finally, repermute the columns */
357  for ( i = 0 ; i < nfree ; i++ )
358  {
359  for ( k = 0 ; k < nfree ; k++ )
360  {
361  hv[per[k]] = matrix2[i][k] ;
362  }
363  for ( k = 0 ; k < nfree ; k++ )
364  {
365  matrix2[i][k] = hv[k] ;
366  }
367  }
368 
369  /* all is well */
370  return 0 ;
371 }
372 
388 static void new_get_mat ( double * xdat,
389  int * xdim,
390  double * ydat,
391  double * wdat,
392  int * ndat,
393  double * fpar,
394  double * epar/*,
395  int * npar */)
396 {
397  double wd ;
398 
399  double yd ;
400  int i, j, n ;
401 
402  for ( j = 0 ; j < nfree ; j++ )
403  {
404  vec[j] = 0.0 ; /* zero sinfo_vector */
405  for ( i = 0 ; i<= j ; i++ ) /* zero sinfo_matrix only on
406  and below diagonal */
407  {
408  matrix1[j][i] = 0.0 ;
409  }
410  }
411  chi2 = 0.0 ; /* reset reduced chi-squared */
412 
413  /* loop through data points */
414  for ( n = 0 ; n < (*ndat) ; n++ )
415  {
416  double wn = wdat[n] ;
417  if ( wn > 0.0 ) /* legal weight ? */
418  {
419  yd=ydat[n] - sinfo_new_gaussian_ellipse(&xdat[(*xdim) * n],fpar) ;
420  sinfo_new_gaussian_ellipse_deriv( &xdat[(*xdim) * n], fpar, epar ) ;
421  chi2 += yd * yd * wn ; /* add to chi-squared */
422  for ( j = 0 ; j < nfree ; j++ )
423  {
424  wd = epar[parptr[j]] * wn ; /* weighted derivative */
425  vec[j] += yd * wd ; /* fill sinfo_vector */
426  for ( i = 0 ; i <= j ; i++ ) /* fill sinfo_matrix */
427  {
428  matrix1[j][i] += epar[parptr[i]] * wd ;
429  }
430  }
431  }
432  }
433 }
434 
435 
463 static int new_get_vec ( double * xdat,
464  int * xdim,
465  double * ydat,
466  double * wdat,
467  int * ndat,
468  double * fpar,
469  double * epar,
470  int * npar )
471 {
472 
473  double dy ;
474  double mii ;
475  double mji ;
476  double mjj ;
477 
478  int i, j, n, r ;
479 
480  /* loop to modify and scale the sinfo_matrix */
481  for ( j = 0 ; j < nfree ; j++ )
482  {
483  mjj = matrix1[j][j] ;
484  if ( mjj <= 0.0 ) /* diagonal element wrong */
485  {
486  return -5 ;
487  }
488  mjj = sqrt( mjj ) ;
489  for ( i = 0 ; i < j ; i++ )
490  {
491  mji = matrix1[j][i] / mjj / sqrt( matrix1[i][i] ) ;
492  matrix2[i][j] = matrix2[j][i] = mji ;
493  }
494  matrix2[j][j] = 1.0 + labda ; /* scaled value on diagonal */
495  }
496 
497  if ( (r = new_inv_mat()) ) /* sinfo_invert sinfo_matrix inlace */
498  {
499  return r ;
500  }
501 
502  for ( i = 0 ; i < (*npar) ; i ++ )
503  {
504  epar[i] = fpar[i] ;
505  }
506 
507  /* loop to calculate correction sinfo_vector */
508  for ( j = 0 ; j < nfree ; j++ )
509  {
510  double dj = 0.0 ;
511  mjj = matrix1[j][j] ;
512  if ( mjj <= 0.0) /* not allowed */
513  {
514  return -7 ;
515  }
516  mjj = sqrt ( mjj ) ;
517  for ( i = 0 ; i < nfree ; i++ )
518  {
519  mii = matrix1[i][i] ;
520  if ( mii <= 0.0 )
521  {
522  return -7 ;
523  }
524  mii = sqrt( mii ) ;
525  dj += vec[i] * matrix2[j][i] / mjj / mii ;
526  }
527  epar[parptr[j]] += dj ; /* new parameters */
528  }
529  chi1 = 0.0 ; /* reset reduced chi-squared */
530 
531  /* loop through the data points */
532  for ( n = 0 ; n < (*ndat) ; n++ )
533  {
534  double wn = wdat[n] ; /* get weight */
535  if ( wn > 0.0 ) /* legal weight */
536  {
537  dy=ydat[n] - sinfo_new_gaussian_ellipse(&xdat[(*xdim) * n],epar);
538  chi1 += wdat[n] * dy * dy ;
539  }
540  }
541  return 0 ;
542 }
543 
544 
545 
594 int sinfo_new_lsqfitd ( double * xdat,
595  int * xdim,
596  double * ydat,
597  double * wdat,
598  int * ndat,
599  double * fpar,
600  double * epar,
601  int * mpar,
602  int * npar,
603  double * tol ,
604  int * its ,
605  double * lab )
606 {
607  int i, n, r ;
608  int itc ; /* fate of fit */
609  int found ; /* fit converged: 1, not yet converged: 0 */
610  int nuse ; /* number of useable data points */
611  double tolerance ; /* accuracy */
612 
613  itc = 0 ; /* fate of fit */
614  found = 0 ; /* reset */
615  nfree = 0 ; /* number of free parameters */
616  nuse = 0 ; /* number of legal data points */
617 
618  if ( *tol < (DBL_EPSILON * 10.0 ) )
619  {
620  tolerance = DBL_EPSILON * 10.0 ; /* default tolerance */
621  }
622  else
623  {
624  tolerance = *tol ; /* tolerance */
625  }
626 
627  labda = fabs( *lab ) * LABFACG ; /* start value for mixing parameter */
628  for ( i = 0 ; i < (*npar) ; i++ )
629  {
630  if ( mpar[i] )
631  {
632  if ( nfree > NPAR ) /* too many free parameters */
633  {
634  return -1 ;
635  }
636  parptr[nfree++] = i ; /* a free parameter */
637  }
638  }
639 
640  if (nfree == 0) /* no free parameters */
641  {
642  return -2 ;
643  }
644 
645  for ( n = 0 ; n < (*ndat) ; n++ )
646  {
647  if ( wdat[n] > 0.0 ) /* legal weight */
648  {
649  nuse ++ ;
650  }
651  }
652 
653  if ( nfree >= nuse )
654  {
655  return -3 ; /* no degrees of freedom */
656  }
657  if ( labda == 0.0 ) /* linear fit */
658  {
659  /* initialize fpar array */
660  for ( i = 0 ; i < nfree ; fpar[parptr[i++]] = 0.0 ) ;
661  new_get_mat ( xdat, xdim, ydat, wdat, ndat, fpar, epar/*, npar */) ;
662  r = new_get_vec ( xdat, xdim, ydat, wdat, ndat, fpar, epar, npar ) ;
663  if ( r ) /* error */
664  {
665  return r ;
666  }
667  for ( i = 0 ; i < (*npar) ; i++ )
668  {
669  fpar[i] = epar[i] ; /* save new parameters */
670  epar[i] = 0.0 ; /* and set errors to zero */
671  }
672  chi1 = sqrt( chi1 / (double) (nuse - nfree) ) ;
673  for ( i = 0 ; i < nfree ; i++ )
674  {
675  if ( (matrix1[i][i] <= 0.0 ) || (matrix2[i][i] <= 0.0) )
676  {
677  return -7 ;
678  }
679  epar[parptr[i]] = chi1 * sqrt( matrix2[i][i] ) /
680  sqrt( matrix1[i][i] ) ;
681  }
682  }
683  else /* non-linear fit */
684  {
685  /*----------------------------------------------------------------
686  * the non-linear fit uses the steepest descent method in combination
687  * with the Taylor method. The mixing of these methods is controlled
688  * by labda. In the outer loop ( called the iteration loop ) we build
689  * the matrix and calculate the correction sinfo_vector. In the
690  inner loop
691  * (called the interpolation loop) we check whether we have obtained a
692  * better solution than the previous one. If so, we leave the inner loop
693  * else we increase lambda ( give more weight to the steepest descent
694  * method) calculate the correction vector and check again. After the
695  * inner loop
696  * we do a final check on the goodness of the fit and if this satisfies
697  * the tolerance we calculate the errors of the fitted parameters.
698  */
699  while ( !found ) /* iteration loop */
700  {
701  if ( itc++ == (*its) ) /* increase iteration counter */
702  {
703  return -4 ;
704  }
705  new_get_mat( xdat, xdim, ydat, wdat, ndat, fpar, epar/*, npar*/ ) ;
706 
707  /*-------------------------------------------------------------
708  * here we decrease labda since we may assume that each iteration
709  * brings us closer to the answer.
710  */
711  if ( labda > LABMING )
712  {
713  labda = labda / LABFACG ; /* decrease labda */
714  }
715  r = new_get_vec ( xdat, xdim, ydat, wdat, ndat, fpar, epar, npar ) ;
716 
717  if ( r ) /* error */
718  {
719  return r ;
720  }
721 
722  while ( chi1 >= chi2 ) /* interpolation loop */
723  {
724  /*-----------------------------------------------------------
725  * The next statement is based on experience, not on the
726  * mathematics of the problem. It is assumed that we have
727  * reached convergence when the pure steepest descent method
728  * does not produce a better solution.
729  */
730  if ( labda > LABMAXG ) /* assume solution found */
731  {
732  break ;
733  }
734  labda = labda * LABFACG ; /* increase mixing parameter */
735  r = new_get_vec(xdat,xdim,ydat,wdat,ndat,fpar,epar,npar) ;
736 
737  if ( r ) /* error */
738  {
739  return r ;
740  }
741  }
742 
743  if ( labda <= LABMAXG ) /* save old parameters */
744  {
745  for ( i = 0 ; i < *npar ; i++ )
746  {
747  fpar[i] = epar[i] ;
748  }
749  }
750  if ( (fabs( chi2 - chi1 ) <= (tolerance * chi1)) ||
751  (labda > LABMAXG) )
752  {
753  /*-----------------------------------------------------------
754  * we have a satisfying solution, so now we need to calculate
755  * the correct errors of the fitted parameters. This we do by
756  * using the pure Taylor
757  * method because we are very close to the real solution.
758  */
759  labda = LABMING ; /* for Taylor solution */
760  new_get_mat(xdat,xdim,ydat,wdat,ndat,fpar,epar/*, npar */) ;
761  r=new_get_vec(xdat,xdim,ydat,wdat,ndat,fpar,epar,npar ) ;
762 
763  if ( r ) /* error */
764  {
765  return r ;
766  }
767  for ( i = 0 ; i < (*npar) ; i++ )
768  {
769  epar[i] = 0.0 ; /* set error to zero */
770  }
771  chi2 = sqrt ( chi2 / (double) (nuse - nfree) ) ;
772 
773  for ( i = 0 ; i < nfree ; i++ )
774  {
775  if ( (matrix1[i][i] <= 0.0) || (matrix2[i][i] <= 0.0) )
776  {
777  return -7 ;
778  }
779  epar[parptr[i]] = chi2 * sqrt( matrix2[i][i] ) /
780  sqrt( matrix1[i][i] ) ;
781  }
782  found = 1 ; /* we found a solution */
783  }
784  }
785  }
786  return itc ; /* return number of iterations */
787 }
788 
820 int
821 sinfo_new_fit_2d_gaussian ( cpl_image * image,
822  double * fit_par,
823  double * derv_par,
824  int * mpar,
825  int lleftx,
826  int llefty,
827  int halfbox_x,
828  int halfbox_y,
829  int* check )
830 {
831  int i, j, n ;
832  int col, row ;
833  int boxi, boxj ;
834  int iters ;
835  int ndata ;
836  int xdim ;
837  int npar ;
838  int its ;
839  double lab ;
840  double tol ;
841  double maxval ;
842  double background ;
843  double amplitude ;
844  float * backarray=NULL ;
845  double M, Mx, My ;
846  double Mxx, Mxy, Myy ;
847  double X0, Y0 ;
848  double xydat[4 *halfbox_x*halfbox_y][XDIMG] ;
849  double zdat[4*halfbox_x*halfbox_y] ;
850  double wdat[4*halfbox_x*halfbox_y] ;
851  double xco, yco ;
852  double value ;
853  double denom ;
854 
855  int llx, lly ;
856  int foundrow ;
857  int foundcol ;
858 
859  int ilx=0;
860  int ily=0;
861  float* pidata=NULL;
862 
863  memset(&wdat[0], 0, (4*halfbox_x*halfbox_y)* sizeof(double));
864 
865  if ( NULL == image )
866  {
867  sinfo_msg_error("no image given") ;
868  return -1 ;
869  }
870  ilx=cpl_image_get_size_x(image);
871  ily=cpl_image_get_size_y(image);
872 
873  if ( NULL == fit_par )
874  {
875  sinfo_msg_error("no fit parameters given") ;
876  return -1 ;
877  }
878  if ( NULL == derv_par )
879  {
880  sinfo_msg_error("no derivatives of fit parameters given") ;
881  return -1 ;
882  }
883  if ( lleftx < 0 || lleftx + 2*halfbox_x >= ilx ||
884  llefty < 0 || llefty + 2*halfbox_y >= ily )
885  {
886  sinfo_msg_error("wrong lower left point of fitting box given!") ;
887  return -1 ;
888  }
889  if ( halfbox_x <= 1 || halfbox_y <= 1 )
890  {
891  sinfo_msg_error("wrong box dimensions given") ;
892  return -1 ;
893  }
894  /* allocate memory */
895  if ( NULL == (backarray = (float*) cpl_calloc(4*halfbox_x+4*halfbox_y,
896  sizeof(float) ) ) )
897  {
898  sinfo_msg_error("could not allocate memory") ;
899  return -1 ;
900  }
901 
902  /* -------------------------------------------------------------------
903  * find the initial estimates for the free parameters
904  */
905 
906  /* first search for the position of the maximum intensity */
907  foundrow = 0 ;
908  foundcol = 0 ;
909  maxval = -SINFO_DBL_MAX ;
910  pidata=cpl_image_get_data_float(image);
911  for ( col = lleftx ; col < lleftx + 2*halfbox_x ; col++ )
912  {
913  for ( row = llefty ; row < llefty + 2*halfbox_y ; row++ )
914  {
915  if ( isnan(pidata[col+row*ilx]) )
916  {
917  continue ;
918  }
919  if ( maxval < pidata[col+row*ilx] )
920  {
921  maxval = pidata[col+row*ilx] ;
922  foundrow = row ;
923  foundcol = col ;
924  }
925  }
926  }
927 
928  if ( foundrow == 0 || foundcol == 0 || maxval <= 0. ||
929  foundrow == ilx-1 || foundcol == ily-1 )
930  {
931  sinfo_msg_warning("no maximum found") ;
932  cpl_free(backarray) ;
933  return -1 ;
934  }
935 
936  /* determine the lower left sinfo_edge of the fitting box, center it
937  on the maximum value */
938  llx = foundcol - halfbox_x ;
939  lly = foundrow - halfbox_y ;
940  if ((foundcol - halfbox_x) > 0) {
941  llx = (foundcol - halfbox_x);
942  } else {
943  llx=1;
944  check++;
945  }
946 
947  if ((foundrow - halfbox_y) > 0) {
948  lly = (foundrow - halfbox_y);
949  } else {
950  lly=1;
951  check++;
952  }
953 
954  if ( ( llx + 2*halfbox_x) < ilx-1 ) {
955  //halfbox_x=halfbox_x;
956  } else {
957  halfbox_x=(int) (ilx-2-llx)/2;
958  check++;
959  }
960 
961  if ( ( lly + 2*halfbox_y) < ily-1 ) {
962  //halfbox_y= halfbox_y;
963  } else {
964  halfbox_y=(int) (ily-2-lly)/2;
965  check++;
966  }
967 
968  if ( llx <= 0 || lly < 0 || llx + 2*halfbox_x >= ilx-1 ||
969  lly + 2*halfbox_y >= ily )
970  {
971  sinfo_msg_error("box does not fit into image") ;
972  cpl_free(backarray) ;
973  return -1 ;
974  }
975 
976  /* determine the zeroth and first order moments of the image
977  within the fitting box */
978  M = Mx = My = 0. ;
979  n = 0 ;
980  boxi = boxj = 0 ;
981  for ( j = lly ; j < lly + 2*halfbox_y ; j++ )
982  {
983  boxj = j - lly ;
984  for ( i = llx ; i < llx + 2*halfbox_x ; i++ )
985  {
986  boxi = i - llx ;
987  if ( !isnan(pidata[i+j*ilx]) )
988  {
989  M += pidata[i+j*ilx] ;
990  Mx += (double)boxi * pidata[i+j*ilx] ;
991  My += (double)boxj * pidata[i+j*ilx] ;
992  /*-----------------------------------------------------------
993  * estimate the amplitude and the background
994  * go through the margins of the fitting box
995  * and calculate the clean mean to
996  * determine the background
997  */
998  if ( i == llx || i == llx + 2*halfbox_x -1 ||
999  j == lly || j == lly + 2*halfbox_y -1 )
1000  {
1001  backarray[n] = pidata[i+j*ilx] ;
1002  n++ ;
1003  }
1004  }
1005  }
1006  }
1007  if ( M <= 0. )
1008  {
1009  sinfo_msg_warning("only negative or zero values") ;
1010  cpl_free(backarray) ;
1011  return -1 ;
1012  }
1013  if ( n < 3 )
1014  {
1015  sinfo_msg_error("not enough data points to calculate background") ;
1016  cpl_free(backarray) ;
1017  return -1 ;
1018  }
1019  /* determine the background as sinfo_median of the surrounding pixels */
1020  if (FLT_MAX==(background=sinfo_new_clean_mean(backarray,n,10.,10.)))
1021  {
1022  sinfo_msg_error("it was not possible to compute the "
1023  "clean mean of the background values") ;
1024  cpl_free(backarray) ;
1025  return -1 ;
1026  }
1027  cpl_free (backarray) ;
1028  /* now calculate the amplitude estimation */
1029  amplitude = maxval - background ;
1030  if ( amplitude < 1e-12 )
1031  {
1032  sinfo_msg_warning("amplitude is too small") ;
1033  return -1 ;
1034  }
1035 
1036  /* determine the center of gravity = centroid */
1037  X0 = Mx / M ;
1038  Y0 = My / M ;
1039  /* if one of the values is outside the fitting box return with error */
1040  if ( X0 <= 0. || Y0 <= 0. || X0 >= 2.*(double)halfbox_x ||
1041  Y0 >= 2.*(double)halfbox_y )
1042  {
1043  sinfo_msg_warning("center of gravity is outside the fitting box!") ;
1044  return -1 ;
1045  }
1046 
1047  /*------------------------------------------------------------------------
1048  * put the data in the 2-d array xydat[][] (pixel position) and zdat[]
1049  * (data values) additionally, determine the second order momentum
1050  */
1051  n = 0 ;
1052  M = Mx = Mxx = My = Myy = Mxy = 0. ;
1053  boxi = boxj = 0 ;
1054  for ( j = lly ; j < lly + 2*halfbox_y ; j++ )
1055  {
1056  boxj = j - lly ;
1057  for ( i = llx ; i < llx + 2*halfbox_x ; i++ )
1058  {
1059  boxi = i - llx ;
1060  value = pidata[i+j*ilx] ;
1061  if ( !isnan(value) )
1062  {
1063  xydat[n][0] = (double) boxi ;
1064  xydat[n][1] = (double) boxj ;
1065  zdat[n] = value ;
1066  wdat[n] = 1. ;
1067  n++ ;
1068 
1069  /* now calculate the moments without background in the
1070  centroid coordinate system */
1071  value -= background ;
1072  xco = (double) boxi - X0 ;
1073  yco = (double) boxj - Y0 ;
1074  M += value ;
1075  Mx += xco * value ;
1076  My += yco * value ;
1077  Mxx += xco * xco * value ;
1078  Myy += yco * yco * value ;
1079  Mxy += xco * yco * value ;
1080  }
1081  }
1082  }
1083  if ( M <= 0. )
1084  {
1085  sinfo_msg_warning("only negative or zero values") ;
1086  return -1 ;
1087  }
1088 
1089  /* ----------------------------------------------------------------
1090  * estimate the fwhm_x and fwhm_y and theta
1091  */
1092 
1093  /* first scale the moments */
1094  /* TODO: why use Mx is later this is never used? */
1095  //Mx /= M ;
1096  //My /= M ;
1097  Mxx /= M ;
1098  Myy /= M ;
1099  Mxy /= M ;
1100 
1101  denom = 2. * (Mxx*Myy - Mxy*Mxy) ;
1102  if ( denom == 0. )
1103  {
1104  sinfo_msg_error("denominator is zero!") ;
1105  return -1 ;
1106  }
1107 
1108  /* now associate the parameter list with the found estimates */
1109  fit_par[0] = X0 ;
1110  fit_par[1] = Y0 ;
1111  fit_par[2] = amplitude ;
1112  fit_par[3] = background ;
1113  fit_par[4] = Myy/denom ;
1114  fit_par[5] = Mxx/denom ;
1115  fit_par[6] = -Mxy/denom ;
1116 
1117  /* convert the moments to ellipse paramters */
1118  if ( 0 > new_gauss2Ellipse (fit_par) )
1119  {
1120  sinfo_msg_warning("gauss2Ellipse does not run!") ;
1121  return -1 ;
1122  }
1123 
1124  /* total number of data points */
1125  ndata = 4 * halfbox_x * halfbox_y ;
1126  xdim = XDIMG ; /* dimension of xydat array */
1127  npar = NPAR ; /* number of parameters in the fit */
1128  its = ITSG ;
1129  lab = LABG ;
1130  tol = TOLG ;
1131  for ( i = 0 ; i < NPAR ; i++ )
1132  {
1133  derv_par[i] = 0. ;
1134  }
1135 
1136  if ( 0 > ( iters = sinfo_new_lsqfitd ( &xydat[0][0],
1137  &xdim,
1138  zdat,
1139  wdat,
1140  &ndata,
1141  fit_par,
1142  derv_par,
1143  mpar,
1144  &npar,
1145  &tol,
1146  &its,
1147  &lab )) )
1148  {
1149  sinfo_msg_warning(" least squares fit failed, error no: %d!", iters) ;
1150  return -1 ;
1151  }
1152 
1153  /* exclude impossible fit results */
1154  if ( fit_par[2] <= 0. || fit_par[4] < 0. || fit_par[5] < 0. )
1155  {
1156  sinfo_msg_error("sorry, some impossible negative fit results!") ;
1157  return -1 ;
1158  }
1159  fit_par[0] += llx ;
1160  fit_par[1] += lly ;
1161  if ( fit_par[0] < llx || fit_par[0] >= llx + 2*halfbox_x ||
1162  fit_par[1] < lly || fit_par[1] >= lly + 2*halfbox_y )
1163  {
1164  sinfo_msg_error("sorry, centroid after the fit "
1165  "outside the fitting box") ;
1166  return -1 ;
1167  }
1168 
1169  /* exchange fwhmx and fwhmy if |theta| is bigger than
1170  pi/4 and subtract pi/2 from theta */
1171  if ( fabs ( fit_par[6] ) > PI_NUMB / 4. )
1172  {
1173  /* first convert angle to smaller than 2 pi */
1174  if ( fabs (fit_par[6]) >= 2. * PI_NUMB )
1175  {
1176  int k = (int) (fit_par[6] / (2.*PI_NUMB)) ;
1177  if ( k > 0 )
1178  {
1179  fit_par[6] -= k*2.*PI_NUMB ;
1180  }
1181  else
1182  {
1183  fit_par[6] += k*2.*PI_NUMB ;
1184  }
1185  }
1186  /* first convert angle to smaller than pi/2 */
1187  if ( fabs (fit_par[6]) > PI_NUMB / 2. )
1188  {
1189  if ( fit_par[6] > 0. )
1190  {
1191  fit_par[6] -= PI_NUMB ;
1192  }
1193  else
1194  {
1195  fit_par[6] += PI_NUMB ;
1196  }
1197  }
1198 
1199  if ( fabs (fit_par[6]) > PI_NUMB / 4. )
1200  {
1201  double temp = fit_par[4] ;
1202  fit_par[4] = fit_par[5] ;
1203  fit_par[5] = temp ;
1204  if ( fit_par[6] < 0. )
1205  {
1206  fit_par[6] += PI_NUMB / 2. ;
1207  }
1208  else
1209  {
1210  fit_par[6] -= PI_NUMB / 2. ;
1211  }
1212  }
1213  }
1214 
1215  return iters ;
1216 }
1217 
1227 cpl_image *
1228 sinfo_crea_gaussian (cpl_image * image,
1229  double * parlist )
1230 {
1231  int col, row ;
1232  cpl_image * retImage ;
1233  double xdat[2] ;
1234  int ilx=0;
1235  int ily=0;
1236  float* podata=NULL;
1237 
1238  if ( image == NULL )
1239  {
1240  sinfo_msg_error("no input image given!") ;
1241  return NULL ;
1242  }
1243  ilx=cpl_image_get_size_x(image);
1244  ily=cpl_image_get_size_y(image);
1245 
1246  if ( parlist == NULL )
1247  {
1248  sinfo_msg_error("no Gaussian parameters given!") ;
1249  return NULL ;
1250  }
1251 
1252  retImage = cpl_image_new (ilx, ily, CPL_TYPE_FLOAT) ;
1253  podata=cpl_image_get_data_float(retImage);
1254  for ( row = 0 ; row < ily ; row++ )
1255  {
1256  for ( col = 0 ; col < ilx ; col++ )
1257  {
1258  xdat[0] = (double) col ;
1259  xdat[1] = (double) row ;
1260  podata[col+row*ilx] = sinfo_new_gaussian_ellipse( xdat , parlist) ;
1261  }
1262  }
1263 
1264  return retImage ;
1265 }
1266 
1274 static int new_gauss2Ellipse ( double * parlist )
1275 {
1276  double a, b, c ;
1277  double ellipseconst ;
1278  double axisX, axisY, phi ;
1279  double p ;
1280 
1281  if ( parlist == NULL )
1282  {
1283  sinfo_msg_error(" no parameters given!\n") ;
1284  return -1 ;
1285  }
1286 
1287  a = parlist[4] ; /* fwhmx */
1288  b = parlist[5] ; /* fwhmy */
1289  c = parlist[6] ; /* theta */
1290 
1291  ellipseconst = 2. * log(2.) ;
1292 
1293  if ( a*b - c*c <= 0. )
1294  {
1295  sinfo_msg_warning("estimates of moments are unusable, "
1296  "they do not make an ellipse!") ;
1297  return -1 ;
1298  }
1299 
1300  if ( a == b )
1301  {
1302  phi = 0. ;
1303  }
1304  else
1305  {
1306  phi = 0.5 * atan( 2. * c / (a-b) ) ;
1307  }
1308 
1309  p = sqrt ( (a-b) * (a-b) + 4. * c*c ) ;
1310 
1311  if ( a > b )
1312  {
1313  axisX = 2. * sqrt ( ellipseconst / (a+b+p) ) ;
1314  axisY = 2. * sqrt ( ellipseconst / (a+b-p) ) ;
1315  }
1316  else
1317  {
1318  axisX = 2. * sqrt ( ellipseconst / (a+b-p) ) ;
1319  axisY = 2. * sqrt ( ellipseconst / (a+b+p) ) ;
1320  }
1321 
1322  parlist[4] = axisX ;
1323  parlist[5] = axisY ;
1324  parlist[6] = phi ;
1325 
1326  return 0 ;
1327 }
1328 
1354 float sinfo_new_determine_conversion_factor ( cpl_imagelist * cube,
1355  float mag,
1356  float exptime,
1357  int llx,
1358  int lly,
1359  int halfbox_x,
1360  int halfbox_y,
1361  int* check )
1362 {
1363  int row, col, i ;
1364  int first_row, first_col ;
1365  int last_row, last_col ;
1366  float factor ;
1367  int mpar[7] ;
1368  double fit_par[7] ;
1369  double derv_par[7] ;
1370 
1371  double sum ;
1372  double xdat[2] ;
1373  cpl_image * summedIm ;
1374 
1375  int ilx=0;
1376  int ily=0;
1377  //int inp=0;
1378 
1379  if ( NULL == cube )
1380  {
1381  sinfo_msg_error(" no cube given!\n") ;
1382  return -FLT_MAX ;
1383  }
1384 
1385  ilx=cpl_image_get_size_x(cpl_imagelist_get(cube,0));
1386  ily=cpl_image_get_size_y(cpl_imagelist_get(cube,0));
1387  //inp=cpl_imagelist_get_size(cube);
1388 
1389  if ( halfbox_x <= 0 || halfbox_y <= 0 ||
1390  2*halfbox_x > ilx || 2*halfbox_y > ily)
1391  {
1392  sinfo_msg_error(" wrong width of halfbox given!") ;
1393  return -FLT_MAX ;
1394  }
1395  if ( exptime <= 0. )
1396  {
1397  sinfo_msg_error(" impossible exposure time given !") ;
1398  return -FLT_MAX ;
1399  }
1400 
1401  /* collapse the cube to be able to do 2D-Gaussian fitting */
1402  if ( NULL == (summedIm = sinfo_new_sum_cube_to_image(cube)) )
1403  {
1404  sinfo_msg_error(" sinfo_averageCubeToImage failed!") ;
1405  return -FLT_MAX ;
1406  }
1407 
1408  /* call the 2D-Gaussian fit routine */
1409  for ( i = 0 ; i < 7 ; i++ )
1410  {
1411  mpar[i] = 1 ;
1412  }
1413  int fitInd ;
1414  if ( -1 == (fitInd = sinfo_new_fit_2d_gaussian(summedIm, fit_par, derv_par,
1415  mpar, llx, lly, halfbox_x,
1416  halfbox_y, check)) )
1417  {
1418  sinfo_msg_warning("sinfo_fit2dGaussian failed!") ;
1419  cpl_image_delete( summedIm) ;
1420  return -FLT_MAX ;
1421  }
1422  cpl_image_delete(summedIm) ;
1423 
1424  /* now integrate the found 2D Gaussian by first
1425  subtracting the background */
1426  if ((fit_par[0] - halfbox_x) < 0) {
1427  first_col=0;
1428  check++;
1429  } else {
1430  first_col=(fit_par[0] - halfbox_x);
1431  }
1432 
1433  if ((fit_par[0] + halfbox_x) < ilx) {
1434  last_col = (fit_par[0] + halfbox_x);
1435  } else {
1436  last_col = (ilx-1) ;
1437  check++;
1438  }
1439 
1440  if ((fit_par[1] - halfbox_y) < 0) {
1441  first_row=0;
1442  check++;
1443  } else {
1444  first_row=(fit_par[1] - halfbox_y) ;
1445  }
1446 
1447  if ((fit_par[1] + halfbox_y) < ily) {
1448  last_row=(fit_par[1] + halfbox_y);
1449  } else {
1450  last_row= (ily-1);
1451  check++;
1452  }
1453 
1454 
1455  if ( first_col < 0 || first_row < 0 || last_col >= ilx || last_row >= ily )
1456  {
1457  sinfo_msg_error("star badly centered in FOV or fitting box too big!") ;
1458  return -FLT_MAX ;
1459  }
1460  sum = 0. ;
1461  for ( row = first_row ; row < last_row ; row++ )
1462  {
1463  for( col = first_col ; col < last_col ; col++ )
1464  {
1465  xdat[0] = (double) col ;
1466  xdat[1] = (double) row ;
1467  sum += (sinfo_new_gaussian_ellipse( xdat, fit_par ) - fit_par[3]) ;
1468  }
1469  }
1470  if ( sum <= 0. )
1471  {
1472  sinfo_msg_error("zero or negative sum of counts!") ;
1473  return -FLT_MAX ;
1474  }
1475  factor = mag / (float)sum * exptime ;
1476  return factor ;
1477 }
1478 
1479 /*--------------------------------------------------------------------------*/
#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