UVES Pipeline Reference Manual  5.5.5b3
uves_physmod_cstacen.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: 2010-09-24 09:32:07 $
23  * $Revision: 1.11 $
24  * $Name: not supported by cvs2svn $
25  */
26 
27 #ifdef HAVE_CONFIG_H
28 # include <config.h>
29 #endif
30 
31 /*---------------------------------------------------------------------------*/
38 /*---------------------------------------------------------------------------*/
39 
40 /* code derived by MIDAS cstacen.c */
41 /*-----------------------------------------------------------------------------
42  Includes
43  ----------------------------------------------------------------------------*/
44 /* definition of the used functions in this module */
45 
46 #include "uves_physmod_cstacen.h"
47 #include <cpl.h>
48 #include <math.h>
49 /*-----------------------------------------------------------------------------
50  Defines
51  ----------------------------------------------------------------------------*/
52 /* Define _POSIX_SOURCE to indicate that this is a POSIX program */
53 /* replaced osmmget by cpl_calloc */
54 /* replaced osmmfree by cpl_free */
55 #define _POSIX_SOURCE 1
56 
57 /* define some macros and constants */
58 
59 #ifndef PI
60 #define PI 3.14159265358979325e0
61 #endif
62 /*
63 #ifndef true
64 #define true 1
65 #define false 0
66 #endif
67 */
68 
69 /* Constants used by the moment centering */
70 
71 #define MINVAL 1.0e-37
72 #define MMXITER 8
73 #define SMALL 1.0e-20
74 
75 
76 /* Constants used by the gaussian centering */
77 
78 #define MAXPAR 4
79 #define IGNORE 2
80 #define NOCONV -1
81 #define OUTSIDE -2
82 #define GMXITER 50
83 #define GCHIMAX 5.0e+16
84 #define GCHIFND 0.005
85 
86 
87 #define MYMIN(a,b) ((a) > (b) ? (b) : (a))
88 #define MYMAX(a,b) ((b) > (a) ? (b) : (a))
89 
92 /*-----------------------------------------------------------------------------
93  Functions prototypes
94  ----------------------------------------------------------------------------*/
95 
96 /*-----------------------------------------------------------------------------
97  Static variables
98  ----------------------------------------------------------------------------*/
99 /*-----------------------------------------------------------------------------
100  Functions code
101  ----------------------------------------------------------------------------*/
102 
103 /*---------------------------------------------------------------------------*/
142 /*---------------------------------------------------------------------------*/
143 
144 int
145 uves_physmod_stacen(float* p_img, int dimx, int dimy, char meth, int* image,
146  float* xout, float* yout, float* xerr, float* yerr,
147  float* xsig, float* ysig, float* xyval, int* stat )
148 
149 {
150 int npix[2], imap[4];
151 float xypos[2], xyerr[2], xysig[2];
152 
153 npix[0] = dimx;
154 npix[1] = dimy;
155 
156 imap[0] = image[0] - 1;
157 imap[1] = image[1] - 1;
158 imap[2] = image[2] - 1;
159 imap[3] = image[3] - 1;
160 
161 /*
162  uves_msg("Input=: npix[0]=%d npix[1]=%d",npix[0],npix[1]);
163  uves_msg("Input=: imap[0]=%d imap[1]=%d imap[2]=%d imap[3]=%d",imap[0],imap[1],imap[2],imap[3]);
164 */
165 
166 *stat = uves_physmod_cstacen(meth, p_img, npix, imap, xypos, xyerr, xysig, xyval );
167 
168 *xout = xypos[0] + 1;
169 *yout = xypos[1] + 1;
170 
171 
172 *xerr = xyerr[0];
173 *yerr = xyerr[1];
174 *xsig = xysig[0];
175 *ysig = xysig[1];
176 /*
177  uves_msg("xout=%f,yout=%f,xerr=%f,yerr=%f,xsig=%f,ysig=%f",
178  *xout,*yout,*xerr,*yerr,*xsig,*ysig);
179  */
180 
181  return 0;
182 }
183 /*---------------------------------------------------------------------------*/
189 static int CGN_NINT(float a){
190  int res=0;
191 
192  res = (a) > 0 ? floor(a+0.5) : -floor(a-0.5);
193  return res;
194 }
195 
196 /*---------------------------------------------------------------------------*/
223 /*---------------------------------------------------------------------------*/
224 
225 
226 static int Ckapsig( float *val, int nval, int iter, float akap,
227  float *cons, float *rms, int *npts )
228 
229 {
230 
231 register int ii=0;
232 register int it=0;
233 register int nr=0;
234 
235 int nr_old=0;
236 
237 float clip=0;
238 float dels=0;
239 float delv=0;
240 float mean=0;
241 float msq=0;
242 float sum=0;
243 float *vsq=NULL;
244 
245 if ( nval < 2 ) return (-1);
246 
247 
248 /* initialize mean value */
249 
250 mean = 0.0;
251 for (ii=0; ii<nval; ii++) mean += val[ii];
252 mean /= (float) (nval);
253 msq = mean * mean;
254 
255 
256 /* initialize RMS */
257 
258 vsq = (float *) cpl_calloc( nval, sizeof( float ));
259 
260 dels = 0.0;
261 for (ii=0; ii<nval; ii++)
262  {
263  vsq[ii] = val[ii] * val[ii];
264  delv = MYMAX( 0.0, vsq[ii] + msq - (2.0 * mean * val[ii]));
265  dels += delv;
266  }
267 
268 *rms = (float) sqrt( MYMAX( MINVAL, dels / (nval-1)));
269 clip = akap * (*rms);
270 
271 
272 /* iterate */
273 
274 nr_old = 0;
275 for (it=0; it<iter; it++)
276  {
277  nr = 0;
278  sum = 0.0;
279  dels = 0.0;
280 
281  for ( ii = 0; ii < nval; ii++ )
282  {
283  if ( fabs( val[ii] - mean ) < (double) clip )
284  {
285  nr++;
286  delv = MYMAX( 0.0, vsq[ii] + msq - 2.0 * mean * val[ii]);
287  dels += delv;
288  sum += val[ii];
289  }
290  }
291 
292  if ( nr <= 2 || nr == nr_old ) goto end_of_it;
293 
294 
295  /* define new rms and mean value */
296 
297  nr_old = nr;
298  *rms = (float) sqrt( MYMAX( MINVAL, dels / (nr-1)));
299  clip = akap * *rms;
300  mean = sum / nr_old;
301  msq = mean * mean;
302  }
303 
304 end_of_it:
305 *cons = mean; /* exit */
306 *npts = nr;
307 
308 cpl_free(vsq);
309 return (0);
310 }
311 
312 
313 /*---------------------------------------------------------------------------*/
328 /*---------------------------------------------------------------------------*/
329 
330 static int MATINV( double (*matrix)[MAXPAR], int nfree )
331 {
332 
333 register int ii=0;
334 register int jj=0;
335 register int kk=0;
336 
337 int evin=0;
338 
339 int per[MAXPAR];
340 double even=0.;
341 double mjk=0.;
342 
343 double hv[MAXPAR];
344 
345 
346 for ( ii = 0; ii < nfree; ii++ ) per[ii] = ii; /* set permutation array */
347 
348 for ( jj = 0; jj < nfree; jj++ ) /* in j-th column, ... */
349  {
350  double rowmax = fabs( matrix[jj][jj] ); /* determine row with ... */
351  int row = jj; /* largest element. */
352  for ( ii = jj + 1; ii < nfree; ii++ )
353  {
354  if ( fabs( matrix[ii][jj] ) > rowmax )
355  {
356  rowmax = fabs( matrix[ii][jj] );
357  row = ii;
358  }
359  }
360 
361  if (fabs(matrix[row][jj]) < SMALL) /* determinant is zero! */
362  return (1);
363 
364  if ( row > jj ) /* if largest element not ...*/
365  {
366  for ( kk = 0; kk < nfree; kk++ ) /* on diagonal, then ... */
367  {
368  even = matrix[jj][kk]; /* permutate rows. */
369  matrix[jj][kk] = matrix[row][kk];
370  matrix[row][kk] = even;
371  }
372  evin = per[jj]; /* keep track of permutation */
373  per[jj] = per[row];
374  per[row] = evin;
375  }
376 
377  even = 1.0 / matrix[jj][jj]; /* modify column */
378  for (ii=0; ii<nfree; ii++)
379  matrix[ii][jj] *= even;
380  matrix[jj][jj] = even;
381  for (kk=0; kk<jj; kk++)
382  {
383  mjk = matrix[jj][kk];
384  for ( ii = 0; ii < jj; ii++ )
385  matrix[ii][kk] -= matrix[ii][jj] * mjk;
386  for ( ii = jj + 1; ii < nfree; ii++ )
387  matrix[ii][kk] -= matrix[ii][jj] * mjk;
388  matrix[jj][kk] = -even * mjk;
389  }
390 
391  for ( kk = jj + 1; kk < nfree; kk++ )
392  {
393  mjk = matrix[jj][kk];
394  for ( ii = 0; ii < jj; ii++ )
395  matrix[ii][kk] -= matrix[ii][jj] * mjk;
396  for ( ii = jj + 1; ii < nfree; ii++ )
397  matrix[ii][kk] -= matrix[ii][jj] * mjk;
398  matrix[jj][kk] = -even * mjk;
399  }
400  }
401 
402 for ( ii = 0; ii < nfree; ii++ ) /* finally, repermute the ...*/
403  {
404  for ( kk = 0; kk < nfree; kk++ ) /* columns. */
405  {
406  hv[per[kk]] = matrix[ii][kk];
407  }
408  for ( kk = 0; kk < nfree; kk++)
409  {
410  matrix[ii][kk] = hv[kk];
411  }
412  }
413 
414 return 0;
415 }
416 
417 /*---------------------------------------------------------------------------*/
428 /*---------------------------------------------------------------------------*/
429 
430 static double ERFCC( double xx )
431 
432 {
433 
434 double z = fabs( xx );
435 double t = 1.0 / (1.0 + (0.5 * z));
436 
437 /*
438 ans = t * exp( -z * z - 1.26551223 + t * ( 1.00002368 + t *
439  ( 0.37409196 + t * ( 0.09678418 + t *
440  (-0.18628806 + t * ( 0.27886807 + t *
441  (-1.13520398 + t * ( 1.48851587 + t *
442  (-0.82215223 + t * 0.17087277 )))))))));
443 
444 
445 the original code above didn't work on Red Hat Linux 5.2
446 with CENTER/GAUSS where the main program is Fortran code
447 neither on Alpha nor on Intel PC (using f2c) ...
448 
449 however it works on SUSE Linux 6.xx on Intel (using g77)
450 
451 therefore this work around which seemed to solve the problem, KB 000522 */
452 
453 double zz = -z * z - 1.26551223 + t * ( 1.00002368 + t *
454  ( 0.37409196 + t * ( 0.09678418 + t *
455  (-0.18628806 + t * ( 0.27886807 + t *
456  (-1.13520398 + t * ( 1.48851587 + t *
457  (-0.82215223 + t * 0.17087277 ))))))));
458 
459 double zozo=0.;
460 if (zz < -500.0)
461  zozo = 0.0;
462 else
463  zozo = exp(zz);
464 
465 double ans = t * zozo;
466 
467 
468 return (xx >= 0.0 ? ans : 2.0 - ans);
469 }
470 
471 /*---------------------------------------------------------------------------*/
481 /*---------------------------------------------------------------------------*/
482 
483 static double GAUSFU( double xx, double *gpar )
484 {
485 
486 double rc=0.;
487 double dd=0.;
488 
489 static int init = TRUE;
490 static double sqrt_2;
491 static double rc1;
492 
493 
494 if ( init )
495  {
496  sqrt_2 = sqrt( 2.0 );
497  rc1 = sqrt(PI)/sqrt_2;
498  init = FALSE;
499  }
500 
501 rc = 1.0 / (sqrt_2 * gpar[2]);
502 dd = xx - gpar[1];
503 dd = ERFCC(rc * (dd - 0.5)) - ERFCC(rc * (dd + 0.5));
504 return ( gpar[3] + rc1 * gpar[0] * gpar[2] * dd );
505 }
506 
507 /*---------------------------------------------------------------------------*/
519 /*---------------------------------------------------------------------------*/
520 
521 static void GAUSDE( double xdat, double *gpar, double *deriv )
522 {
523 double temp=0.;
524 double tempp=0.;
525 double x1=0.;
526 double x2=0.;
527 double zz=0.;
528 
529 
530 
531 static double sqrt_2;
532 
533 register int jj=0;
534 static int init = TRUE;
535 
536 
537 if ( init )
538  {
539  sqrt_2 = sqrt( 2.0 );
540  init = FALSE;
541  }
542 
543 temp = sqrt_2 * gpar[2];
544 tempp = xdat - gpar[1];
545 x1 = (tempp - 0.5) / temp;
546 x2 = (tempp + 0.5) / temp;
547 zz = tempp / gpar[2] ;
548 
549 if ( ((zz * zz) - 50.0) < 0.0 )
550  {
551  deriv[0] = (GAUSFU( xdat, gpar ) - gpar[3]) / gpar[0];
552  double zx = (-x1) * x1;
553  double dv1=0.;
554  if ( zx < -200.0) /* zx always < 0 */
555  dv1 = 0.0; /* e**(-200) is = 0.0 ... */
556  else
557  dv1 = exp(zx);
558  zx = (-x2) * x2;
559  double dv2=0.;
560  if ( zx < -200.0)
561  dv2 = dv1;
562  else
563  dv2 = dv1 - exp(zx);
564  deriv[1] = gpar[0] * dv2;
565 
566  /* for (x1 * x1) > 400 we got floating point exceptions on DEC Alpha
567  deriv[1] = gpar[0] * (exp( -x1 * x1 ) - exp( -x2 * x2 ));
568  */
569 
570  deriv[2] = deriv[1] * zz;
571  }
572 else
573  for (jj=0; jj<3; jj++) deriv[jj] = 0.0;
574 
575 deriv[3] = 1.0;
576 }
577 
578 /*---------------------------------------------------------------------------*/
594 /*---------------------------------------------------------------------------*/
595 
596 static float FCHIS(double *data,int ndim,int nfree,int mode,double *dfit) {
597 
598 register int ii=0;
599 
600 
601 
602 
603 
604 if ( nfree > 0 )
605  {
606  double chisq = 0.0;
607 
608  for (ii=0; ii<ndim; ii++)
609  {
610  double weight=0.;
611  if ( mode < 0 )
612  {
613  if ( *data < 0 )
614  weight = -1. / *data;
615  else if ( *data == 0 )
616  weight = 1.0;
617  else
618  weight = 1. / *data;
619  }
620  else
621  weight = 1.0;
622 
623  double diff = (*data) - (*dfit);
624  data++; dfit++;
625  chisq += weight * diff * diff;
626  }
627  return (chisq / nfree);
628  }
629 
630 else
631  return 0.0;
632 }
633 
634 /*---------------------------------------------------------------------------*/
662 static int LSQFIT( double *xdat, double *data, int ndim,
663  double *gpar, float *lamda, double *dfit,
664  double *chisqr, double *sigma )
665 
666 {
667 
668 register int icnt=0;
669 register int ii=0;
670 register int jj=0;
671 register int kk=0;
672 
673 int nfree=0;
674 
675 double chisq1=0;
676 double b[MAXPAR], beta[MAXPAR], deriv[MAXPAR],
677  array[MAXPAR][MAXPAR], alpha[MAXPAR][MAXPAR];
678 
679 
680 nfree = ndim - MAXPAR;
681 *sigma = 0.0;
682 if ( nfree < 1 || fabs( (double) *gpar ) < SMALL ) return (1);
683 
684 
685 /* evaluate ALPHA and BETA matrices */
686 
687 for (ii=0; ii<MAXPAR; ii++)
688  {
689  beta[ii] = 0.0;
690  for (jj=0; jj<=ii; jj++) alpha[ii][jj] = 0.0;
691  }
692 
693 for (ii=0; ii<ndim; ii++)
694  {
695  GAUSDE( xdat[ii], gpar, deriv ); /* here we divide by gpar[1] */
696 
697  for (jj=0; jj<MAXPAR; jj++)
698  {
699  beta[jj] += (data[ii] - GAUSFU( xdat[ii], gpar )) * deriv[jj];
700  for (kk=0; kk<=jj; kk++)
701  alpha[jj][kk] += deriv[jj] * deriv[kk];
702  }
703  }
704 
705 for (ii=0; ii<MAXPAR; ii++)
706  {
707  for (jj=0; jj<=ii; jj++)
708  alpha[jj][ii] = alpha[ii][jj];
709  }
710 
711 
712 /* invert matrix */
713 
714 if ( *lamda < SMALL)
715  {
716  if (MATINV(alpha,MAXPAR) == 1) return (2); /* determinant -> 0.0 */
717 
718  *sigma = MYMAX( 0.0, alpha[1][1] );
719  }
720 
721 else /* evaluate chi square at starting point */
722  {
723  for (ii=0; ii<ndim; ii++)
724  dfit[ii] = GAUSFU( xdat[ii], gpar );
725 
726  chisq1 = FCHIS( data, ndim, nfree, 0, dfit );
727 
728  icnt = 0; /* invert matrix */
729 loop:
730  for ( jj = 0; jj < MAXPAR; jj++ )
731  {
732  for ( kk = 0; kk < MAXPAR; kk++ )
733  {
734  if (fabs( alpha[jj][jj] ) < 1.e-15 || fabs( alpha[kk][kk] ) < 1.e-15)
735  return 2;
736  array[jj][kk] = alpha[jj][kk] /
737  sqrt( alpha[jj][jj] * alpha[kk][kk] ) ;
738  }
739  array[jj][jj] = 1.0 + *lamda;
740  }
741 
742  (void) MATINV( array, MAXPAR );
743 
744  for ( jj = 0; jj < MAXPAR; jj++ )
745  {
746  b[jj] = gpar[jj] ;
747  for ( kk = 0; kk < MAXPAR ; kk++ )
748  {
749  b[jj] += beta[kk] * array[jj][kk] /
750  sqrt( alpha[jj][jj] * alpha[kk][kk] );
751  }
752  }
753 
754 /* if chi square increased, increase LAMDA and try again */
755 
756  for (ii=0; ii<ndim; ii++)
757  dfit[ii] = GAUSFU( xdat[ii], b );
758 
759  *chisqr = FCHIS( data, ndim, nfree, 0, dfit );
760 
761  if ( chisq1 - *chisqr < 0.0 )
762  {
763  if (++icnt < 60)
764  {
765  *lamda *= 10;
766  goto loop;
767  }
768  else
769  return (2);
770  }
771 
772  for (jj=0; jj<MAXPAR; jj++) gpar[jj] = b[jj];
773  *lamda /= 10.0;
774  }
775 
776 return 0;
777 }
778 
779 
780 /*---------------------------------------------------------------------------*/
800 static void Crhox( float *p_img, int *npix, int *image,
801  int *lnew, double *krx )
802 
803 {
804 register int nxdim=0;
805 register int ix=0;
806 register int iy=0;
807 
808 int nrx=0;
809 int nry=0;
810 
811 
812 
813 
814 
815 
816 /* original FORTRAN code:
817 
818  IXA = IMAP(1) IMAP(1-4) => image[0-3]
819  IXE = IMAP(2)
820  IYA = IMAP(3)
821  JYA = IYA + JY - 1 JY => lnew[0]
822  JYE = IYA + LY - 1 LY => lnew[1]
823  M = 1
824 
825  DO 200 J=IXA,IXE
826  ISUM = 0.D0
827  DO 100 K=JYA,JYE
828  ISUM = ISUM + AIMG(J,K)
829 100 CONTINUE
830  KRX(M) = ISUM
831  M = M + 1
832 200 CONTINUE
833 */
834 
835 
836 nrx = image[1] - image[0] + 1;
837 nry = lnew[1] - lnew[0] + 1;
838 nxdim = *npix;
839 p_img += nxdim * (image[2] + lnew[0]);
840 
841 for (ix=0; ix<nrx; ix++)
842  {
843  double sum = 0.0;
844  for (iy=0; iy<nry*nxdim; iy+=nxdim) sum += p_img[iy];
845  p_img++;
846  *krx++ = sum;
847  }
848 }
849 
850 
851 
852 /*---------------------------------------------------------------------------*/
872 static void Crhoy( float *p_img, int *npix, int *image,
873  int *lnew, double *kry )
874 
875 {
876 register int nxdim=0;
877 register int ix=0;
878 register int iy=0;
879 
880 int nrx=0;
881 int nry=0;
882 
883 
884 
885 
886 
887 /* original FORTRAN code:
888 
889  IXA = IMAP(1) IMAP(1-4) => image[0-3]
890  IYA = IMAP(3)
891  IYE = IMAP(4)
892  JXA = IXA + JX - 1 JX => lnew[0]
893  JXE = IXA + LX - 1 LX => lnew[1]
894 
895  M = 1
896  DO 200 J=IYA,IYE
897  ISUM = 0.D0
898  DO 100 K=JXA,JXE
899  ISUM = ISUM + AIMG(K,J)
900 100 CONTINUE
901  KRY(M) = ISUM
902  M = M + 1
903 200 CONTINUE
904 
905 */
906 
907 nrx = lnew[1] - lnew[0] + 1;
908 nry = image[3] - image[2] + 1;
909 nxdim = *npix;
910 p_img += (nxdim * image[2]) + (image[0] + lnew[0]);
911 
912 for (iy=0; iy<nry; iy++)
913  {
914  double sum = 0.0;
915  for (ix=0; ix<nrx; ix++) sum += p_img[ix];
916  p_img += nxdim;
917  *kry++ = sum;
918  }
919 }
920 
921 /*---------------------------------------------------------------------------*/
945 static int Cserch( double *marg, int ndim, int ign,
946  int *lmin, int *lmax, float *s_cent, float *s_width )
947 
948 {
949 
950 register int ii=0;
951 
952 int ql=0;
953 int ibgn=0;
954 int icrowd=0;
955 int iend=0;
956 int imax=0;
957 int imin=0;
958 int indx=0;
959 
960 double dxk=0.;
961 double diff=0.;
962 
963 
964 double sum=0.;
965 double *work=NULL;
966 
967 
968 ibgn = ign;
969 iend = ndim - ign -1; /* ojo */
970 
971 /* create workspace to store the derivative of the marginal data */
972 
973 work = (double *) cpl_calloc( ndim , sizeof( double ));
974 
975 
976 /* find maximum and minimum derivative of MARG */
977 
978 imin = imax = 0;
979 double drmn = 0., drmx = 0.0;
980 for (ii = ibgn; ii < iend; ii++ )
981  {
982  diff = marg[ii+1] - marg[ii-1];
983  work[ii] = marg[ii+2] - marg[ii-2] + (2 * diff);
984  if ( work[ii] >= drmx )
985  {
986  drmx = work[ii];
987  imax = ii;
988  }
989  if (
990  work[ii] <= drmn )
991  { drmn = work[ii];
992  imin = ii;
993  }
994  }
995 
996 
997 /* crowded ? */
998 
999 icrowd = 0;
1000 if (imin <= imax ) /* bright source to the left, compute right a new minima */
1001  {
1002  if ( ndim - imax > imin )
1003  {
1004  icrowd = -1;
1005  drmn = drmx;
1006  for ( ii = imax+1; ii < iend; ii++ )
1007  {
1008  if ( work[ii] < drmn )
1009  {
1010  drmn = work[ii];
1011  imin = ii;
1012  }
1013  }
1014  }
1015  else /* bright source to the right, compute left a new maxima */
1016  {
1017  icrowd = 1;
1018  drmx = drmn;
1019  for ( ii = ibgn; ii < imin; ii++ )
1020  {
1021  if ( work[ii] >= drmx )
1022  {
1023  drmx = work[ii];
1024  imax = ii;
1025  }
1026  }
1027  }
1028  }
1029 
1030 
1031 /* compute estimates of image centre and width */
1032 
1033 *s_cent = ((float)(imax + imin)) * 0.5;
1034 *s_width = imin - imax;
1035 
1036 sum = 0.0;
1037 for ( ii = imax; ii <= imin; ii++ ) sum += work[ii];
1038 
1039 diff = drmx - drmn;
1040 if ( fabs(diff) > SMALL)
1041  {
1042  dxk = sum * *s_width / ( (*s_width+1.0)*diff );
1043  *s_cent += dxk;
1044  }
1045 *s_width /= 2;
1046 indx = CGN_NINT(*s_cent);
1047 if (indx < 0)
1048  {
1049  *s_cent = 0.0;
1050  indx = 0;
1051  }
1052 else if (indx >= ndim)
1053  {
1054  *s_cent = (float)(ndim-1);
1055  indx = ndim-1;
1056  }
1057 
1058 
1059 /* find low- (left-) side local minimum */
1060 
1061 ql = indx - 2;
1062 *lmin = 0;
1063 if (ql < 2) goto next_step;
1064 
1065 low_loop:
1066 ql --;
1067 if (ql <= 0 ) goto next_step;
1068 if (ql == 1) goto lo5;
1069 if (ql == 2) goto lo4;
1070 if (ql == 3) goto lo3;
1071 
1072 if (marg[ql] > marg[ql-4]) goto low_loop;
1073 lo3:
1074 if (marg[ql] > marg[ql-3]) goto low_loop;
1075 lo4:
1076 if (marg[ql] > marg[ql-2]) goto low_loop;
1077 lo5:
1078 if (marg[ql] > marg[ql-1]) goto low_loop;
1079 
1080 *lmin = ql + 1;
1081 
1082 
1083 /* find high- (right-) side local minimum */
1084 
1085 next_step:
1086 ql = indx + 2;
1087 *lmax = ndim - 1;
1088 ii = ndim - ql;
1089 if (ii < 3) goto end_of_it;
1090 
1091 hi_loop:
1092 ql ++ ;
1093 ii = ndim - ql;
1094 if (ii == 1 ) goto end_of_it;
1095 if (ii == 2) goto hi5;
1096 if (ii == 3) goto hi4;
1097 if (ii == 4) goto hi3;
1098 
1099 if (marg[ql] > marg[ql+4]) goto hi_loop;
1100 hi3:
1101 if (marg[ql] > marg[ql+3]) goto hi_loop;
1102 hi4:
1103 if (marg[ql] > marg[ql+2]) goto hi_loop;
1104 hi5:
1105 if (marg[ql] > marg[ql+1]) goto hi_loop;
1106 *lmax = ql - 1;
1107 
1108 end_of_it:
1109 (void) cpl_free( (char *) work );
1110 return (icrowd);
1111 
1112 }
1113 
1114 /*---------------------------------------------------------------------------*/
1152 /*---------------------------------------------------------------------------*/
1153 
1154 /* ------------------------------------------*/
1155 /* here starts the code of the main function */
1156 /* ------------------------------------------*/
1157 
1158 
1159 int uves_physmod_cstacen(char meth, float* p_img, int* npix, int* image,
1160  float* xypos, float* xyerr, float* xysig, float* xyval )
1161 
1162 {
1163 
1164 register int it=0;
1165 register int ix=0;
1166 register int iy=0;
1167 
1168 int bgnr=0;
1169 int indx=0;
1170 int indy=0;
1171 int istat=0;
1172 int nrx=0;
1173 int nry=0;
1174 int nval=0;
1175 int ifram[4];
1176 
1177 float bgval=0.;
1178 float clip=0.;
1179 float rms=0.;
1180 float xmom=0.;
1181 float ymom=0.;
1182 float source=0.;
1183 float sumi=0.;
1184 float xold=0.;
1185 float yold=0.;
1186 
1187 float *p_bgn=NULL;
1188 float *p_edge=NULL;
1189 
1190 
1191 istat = 0; /* initialize */
1192 p_bgn = p_img;
1193 for (ix=0; ix<4; ix++)
1194  ifram[ix] = image[ix] + 1; /* 1 ---> ndim */
1195 nrx = ifram[1] - ifram[0] + 1;
1196 nry = ifram[3] - ifram[2] + 1;
1197 
1198 xypos[0] = (ifram[0] + ifram[1]) * 0.5; /* init to center pixel */
1199 xypos[1] = (ifram[2] + ifram[3]) * 0.5;
1200 xyerr[0] = xyerr[1] = 0.0;
1201 xysig[0] = xysig[1] = 0.0;
1202 *xyval = 0.0;
1203 
1204 /* MOMENT centering */
1205 
1206 if ( meth != 'G' && meth != 'g' )
1207  {
1208  int kk, istr, iend;
1209 
1210  xold = yold = -1.0; /* find bgval and rms from edge pixels */
1211 
1212  p_img += (ifram[0] - 1) + (npix[0] * (ifram[2] - 1));
1213 
1214  /* collect edge pixels */
1215 
1216  if (nry > 1)
1217  {
1218  nval = (2 * nrx) + (2 * (nry-2));
1219  p_edge = (float *) cpl_calloc( nval , sizeof( float ));
1220 
1221  for (ix=0; ix<nrx;ix++)
1222  *p_edge++ = p_img[ix];
1223 
1224  p_img += *npix;
1225  for (iy=0; iy<(nry-2); iy++)
1226  {
1227  *p_edge++ = p_img[0];
1228  *p_edge++ = p_img[nrx - 1];
1229  p_img += npix[0];
1230  }
1231  for (ix=0; ix<nrx; ix++) *p_edge++ = p_img[ix];
1232  }
1233  else
1234  {
1235  nval = nrx;
1236  p_edge = (float *) cpl_calloc( nval , sizeof( float ));
1237 
1238  for (ix=0; ix<nrx;ix++)
1239  *p_edge++ = p_img[ix];
1240  }
1241 
1242  p_img = p_bgn;
1243 
1244  p_edge -= nval;
1245  (void) Ckapsig( p_edge, nval, 5, 2.0, &bgval, &rms, &bgnr );
1246  (void) cpl_free( (char *) p_edge );
1247 
1248 
1249  /* calculate moment for pixel values > 3 * RMS above BGVAL */
1250 
1251  clip = 3.0 * rms;
1252 
1253  for (it=0; it<MMXITER; it++) /* iteration loop */
1254  {
1255  sumi = xmom = ymom = 0.0;
1256  p_img += ifram[0] - 1 + (npix[0] * (ifram[2] - 1));
1257  for (iy=0; iy<nry; iy++)
1258  {
1259  for (ix=0; ix<nrx; ix++)
1260  {
1261  if ( (source = p_img[ix] - bgval) > clip )
1262  {
1263  sumi += source;
1264  xmom += source * (ifram[0] + ix);
1265  ymom += source * (ifram[2] + iy);
1266  }
1267  }
1268  p_img += npix[0];
1269  }
1270  p_img = p_bgn; /* reset to start of array */
1271 
1272  if ((nrx < 3) || (nry < 3))
1273  {
1274  xysig[0] = nrx;
1275  xysig[1] = nry;
1276  istat = 1;
1277  if ( sumi > 0.0 )
1278  {
1279  xypos[0] = xmom / sumi;
1280  xypos[1] = ymom / sumi;
1281  }
1282  else
1283  {
1284  istat = 2;
1285  xypos[0] = (ifram[0] + ifram[1]) * 0.5;
1286  xypos[1] = (ifram[2] + ifram[3]) * 0.5;
1287  }
1288  indx = CGN_NINT(xypos[0]-1);
1289  indy = CGN_NINT(xypos[1]-1);
1290  *xyval = p_img[indx + ((*npix) * indy)];
1291  goto end_of_iter; /* EXIT iteration loop */
1292  }
1293 
1294  if (sumi > 0.0) /* only positive sources */
1295  {
1296  xypos[0] = xmom / sumi;
1297  xypos[1] = ymom / sumi;
1298  xysig[0] = nrx;
1299  xysig[1] = nry;
1300 
1301  if ( xold == xypos[0] && yold == xypos[1] )
1302  {
1303  int nr = 0;
1304  double xdif, ydif, xrms, yrms;
1305 
1306  xrms = yrms = sumi = 0.0;
1307  p_img += ifram[0] - 1 + (npix[0] * (ifram[2] - 1));
1308  for (iy=0; iy<nry; iy++ )
1309  {
1310  for (ix=0; ix<nrx; ix++)
1311  {
1312  if ( (source = p_img[ix] - bgval) > clip )
1313  {
1314  xdif = (ifram[0] + ix) - xypos[0];
1315  ydif = (ifram[2] + iy) - xypos[1];
1316  xrms += fabs( source * xdif *xdif );
1317  yrms += fabs( source * ydif *ydif );
1318  sumi += source;
1319  nr++;
1320  }
1321  }
1322  p_img += npix[0];
1323  }
1324  p_img = p_bgn;
1325 
1326  indx = CGN_NINT(xypos[0]-1) + (npix[0] * CGN_NINT(xypos[1]-1));
1327  *xyval = p_img[indx];
1328  xysig[0] = (float) sqrt(xrms /(sumi+ *xyval - bgval));
1329  xysig[1] = (float) sqrt(yrms /(sumi+ *xyval - bgval));
1330  xyerr[0] = (float) (xysig[0] / sqrt( (double) (nr - 1)));
1331  xyerr[1] = (float) (xysig[1] / sqrt( (double) (nr - 1)));
1332  goto end_of_iter; /* succesful return */
1333  }
1334 
1335 
1336  xold = xypos[0];
1337  yold = xypos[1];
1338  }
1339  else
1340  {
1341  istat = 2;
1342  xypos[0] = (ifram[0] + ifram[1]) * 0.5;
1343  xypos[1] = (ifram[2] + ifram[3]) * 0.5;
1344  indx = CGN_NINT(xypos[0]-1);
1345  indy = CGN_NINT(xypos[1]-1);
1346  *xyval = p_img[indx + ((*npix) * indy)];
1347  goto end_of_iter; /* EXIT iteration loop */
1348  }
1349 
1350 
1351  /* crowded or weak source conditions */
1352 
1353  indx = CGN_NINT(xypos[0]-1) + (npix[0] * CGN_NINT(xypos[1]-1));
1354  if ( (*xyval = p_img[indx] - bgval) <= clip )
1355  {
1356  xysig[0] = xysig[1] = 0.0;
1357  istat = 1;
1358  goto end_of_iter; /* EXIT iteration loop */
1359  }
1360 
1361 
1362  /* find extent of source i.e. delete spikes, etc. */
1363 
1364  ix = CGN_NINT( xypos[0] ); /* ix, iy = 1,2,... */
1365  iy = CGN_NINT( xypos[1] );
1366  kk = npix[0] * (iy - 1);
1367  istr = ifram[0];
1368  source = p_img[ix-1 + kk] - bgval;
1369  while ( source > clip && ix >= istr )
1370  {
1371  ifram[0] = ix;
1372  source = p_img[ix-1 + kk] - bgval;
1373  ix --;
1374  }
1375 
1376  ix = CGN_NINT( xypos[0] );
1377  iend = ifram[1];
1378  source = p_img[ix-1 + kk] - bgval;
1379  while ( source > clip && ix <= iend )
1380  {
1381  ifram[1] = ix;
1382  source = p_img[ix-1 + kk] -bgval;
1383  ix ++;
1384  }
1385 
1386  ix = CGN_NINT( xypos[0] );
1387  istr = ifram[2];
1388  source = p_img[ix-1 + kk] - bgval;
1389  while ( source > clip && iy >= istr )
1390  {
1391  ifram[2] = iy;
1392  source = p_img[ix-1 + (*npix *(iy-1))] -bgval;
1393  iy --;
1394  }
1395 
1396  iy = CGN_NINT( xypos[1] );
1397  iend = ifram[3];
1398  source = p_img[ix-1 + kk] - bgval;
1399  while ( source > clip && iy <= iend )
1400  {
1401  ifram[3] = iy;
1402  source = p_img[ix-1 + (*npix *(iy-1))] -bgval;
1403  iy++;
1404  }
1405  nrx = ifram[1] - ifram[0] + 1;
1406  nry = ifram[3] - ifram[2] + 1;
1407  }
1408 
1409  istat = 3; /* iteration failed */
1410 
1411 end_of_iter:
1412  xypos[0] --;
1413  xypos[1] --;
1414  }
1415 
1416 
1417 /* GAUSSIAN centering */
1418 
1419 else
1420  {
1421  register int ii;
1422  int found, ierr, xlim[2], ylim[2], lnew[2];
1423  float lamda, xcent, ycent, xwidth, ywidth;
1424  double chisqr, oldchi, sigma, *krx, *kry, *gfit, *xpos, *yfit,
1425  gpar[MAXPAR];
1426 
1427 
1428  /* construct two marginal distibutions (in pixel coordinates!) */
1429 
1430  lnew[0] = (nry / 4); /* in C notation, from 0 ... */
1431  lnew[1] = nry - (nry / 4) - 1;
1432 
1433 
1434 /* Take care of 1-dim case */
1435 
1436  if (nry == 1)
1437  {
1438  krx = (double *) cpl_calloc(nrx , sizeof(double));
1439  Crhox(p_img,npix,image,lnew,krx);
1440  ierr = Cserch(krx,nrx,IGNORE,xlim,xlim+1,&xcent,&xwidth);
1441 
1442  /* store the data of the fit */
1443 
1444  nval = xlim[1] - xlim[0] + 1;
1445  xpos = (double *) cpl_calloc( nval , sizeof( double ));
1446  yfit = (double *) cpl_calloc( nval , sizeof( double ));
1447  gfit = (double *) cpl_calloc( nval , sizeof( double ));
1448  for (ii=0; ii<nval; ii++)
1449  {
1450  xpos[ii] = xlim[0] + ii;
1451  yfit[ii] = krx[xlim[0] + ii];
1452  }
1453 
1454  /* set parameters for LSQFIT (old FITINTE) */
1455 
1456  lamda = 0.001;
1457  chisqr = GCHIMAX;
1458  gpar[0] = krx[CGN_NINT(xcent)];
1459  gpar[1] = xcent;
1460  gpar[2] = xwidth;
1461  gpar[3] = (krx[xlim[0]] + krx[xlim[1]]) / 2;
1462  (void) cpl_free( (char *) krx );
1463 
1464  it = 0;
1465  found = FALSE;
1466  while ( ! found && it++ < GMXITER )
1467  {
1468  oldchi = chisqr;
1469  ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma);
1470  if ( ierr != 0 )
1471  {
1472  found = NOCONV;
1473  istat = 3;
1474  }
1475  else if ( (oldchi - chisqr)/ chisqr < GCHIFND )
1476  found = TRUE;
1477  }
1478 
1479  /* Is the source still in the image and has it a resonable shape? */
1480 
1481  lamda = 0.0;
1482  ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma);
1483  if ( ierr != 0 )
1484  {
1485  found = NOCONV;
1486  istat = 3;
1487  }
1488  else
1489  {
1490  sumi = (float)(gpar[1] + image[0]);
1491  indx = CGN_NINT(sumi);
1492  if ( indx < 0 || indx >= *npix )
1493  {
1494  found = OUTSIDE;
1495  istat = 2;
1496  }
1497  }
1498  (void) cpl_free( (char *) xpos );
1499  (void) cpl_free( (char *) yfit );
1500  (void) cpl_free( (char *) gfit );
1501 
1502  if ( found == TRUE )
1503  {
1504  xypos[0] = sumi;
1505  xypos[1] = 0;
1506  xysig[0] = (float) gpar[2];
1507  xyerr[0] = (float) sqrt( sigma * chisqr );
1508  indx = CGN_NINT( xypos[0]);
1509  *xyval = p_img[indx];
1510  }
1511  }
1512 
1513 
1514 /* Take care of 2-dim case */
1515 
1516  else
1517  {
1518  krx = (double *) cpl_calloc( nrx , sizeof( double ));
1519  kry = (double *) cpl_calloc( nry , sizeof( double ));
1520 
1521  /* Compute and search X-marginal & Y-marginal */
1522 
1523  Crhox( p_img, npix, image, lnew, krx );
1524  ierr = Cserch( krx, nrx, IGNORE, xlim, xlim+1, &xcent, &xwidth );
1525  lnew[0] = MYMAX( xlim[0], CGN_NINT(xcent - (2 * xwidth)));
1526  lnew[1] = MYMIN( xlim[1], CGN_NINT(xcent + (2 * xwidth)));
1527 
1528  Crhoy( p_img, npix, image, lnew, kry );
1529  ierr = Cserch( kry, nry, IGNORE, ylim, ylim+1, &ycent, &ywidth );
1530  lnew[0] = MYMAX( ylim[0], CGN_NINT(ycent - (2 * ywidth)));
1531  lnew[1] = MYMIN( ylim[1], CGN_NINT(ycent + (2 * ywidth)));
1532 
1533  Crhox( p_img, npix, image, lnew, krx );
1534  ierr = Cserch( krx, nrx, IGNORE, xlim, xlim+1, &xcent, &xwidth );
1535  lnew[0] = MYMAX( xlim[0], CGN_NINT(xcent - (2 * xwidth)));
1536  lnew[1] = MYMIN( xlim[1], CGN_NINT(xcent + (2 * xwidth)));
1537 
1538  Crhoy( p_img, npix, image, lnew, kry );
1539  ierr = Cserch( kry, nry, IGNORE, ylim, ylim+1, &ycent, &ywidth );
1540 
1541  /* fit a gaussian to the source along the X-axis */
1542 
1543  nval = xlim[1] - xlim[0] + 1;
1544  xpos = (double *) cpl_calloc( nval , sizeof( double ));
1545  yfit = (double *) cpl_calloc( nval , sizeof( double ));
1546  gfit = (double *) cpl_calloc( nval , sizeof( double ));
1547  for (ii=0; ii<nval; ii++)
1548  {
1549  xpos[ii] = xlim[0] + ii;
1550  yfit[ii] = krx[xlim[0] + ii];
1551  }
1552 
1553  /* set parameters for LSQFIT */
1554 
1555  lamda = 0.001;
1556  chisqr = GCHIMAX;
1557  gpar[0] = krx[CGN_NINT( xcent )];
1558  gpar[1] = xcent;
1559  gpar[2] = xwidth;
1560  gpar[3] = (krx[xlim[0]] + krx[xlim[1]]) / 2;
1561  (void) cpl_free( ( char *) krx );
1562 
1563  it = 0;
1564  found = FALSE;
1565  while ( ! found && it++ < GMXITER )
1566  {
1567  oldchi = chisqr;
1568  ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma);
1569  if ( ierr != 0 || gpar[2] <= 0.0 )
1570  {
1571  found = NOCONV;
1572  istat = 3;
1573  }
1574  else if ( (oldchi - chisqr)/ chisqr < GCHIFND )
1575  found = TRUE;
1576  }
1577 
1578  /* Is the source still in the image and has it a resonable shape? */
1579 
1580  lamda = 0.0;
1581  ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma);
1582  if ( ierr != 0 )
1583  {
1584  found = NOCONV;
1585  istat = 3;
1586  }
1587  else
1588  {
1589  sumi = (float)(gpar[1] + image[0]);
1590  indx = CGN_NINT(sumi);
1591  if ( indx < 0 || indx >= *npix )
1592  {
1593  found = OUTSIDE;
1594  istat = 2;
1595  }
1596  }
1597 
1598  (void) cpl_free( (char *) xpos );
1599  (void) cpl_free( (char *) yfit );
1600  (void) cpl_free( (char *) gfit );
1601 
1602  if ( found == TRUE )
1603  {
1604  xypos[0] = sumi;
1605  xysig[0] = (float) gpar[2];
1606  xyerr[0] = (float) sqrt( sigma * chisqr );
1607 
1608  /* x-dir o.k. - now fit a gaussian to the source along the Y-axis */
1609 
1610  nval = ylim[1] - ylim[0] + 1;
1611  xpos = (double *) cpl_calloc( nval , sizeof( double ));
1612  yfit = (double *) cpl_calloc( nval , sizeof( double ));
1613  gfit = (double *) cpl_calloc( nval , sizeof( double ));
1614 
1615  for (ii=0; ii<nval; ii++)
1616  {
1617  xpos[ii] = ylim[0] + ii;
1618  yfit[ii] = kry[ylim[0] + ii];
1619  }
1620 
1621  /* set parameters for LSQFIT */
1622 
1623  lamda = 0.001;
1624  chisqr = GCHIMAX;
1625  gpar[0] = kry[CGN_NINT( ycent )];
1626  gpar[1] = ycent;
1627  gpar[2] = ywidth;
1628  gpar[3] = (kry[ylim[0]] + kry[ylim[1]]) / 2;
1629 
1630  it = 0;
1631  found = FALSE;
1632  while ( ! found && it++ < GMXITER )
1633  {
1634  oldchi = chisqr;
1635  ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma );
1636  if ( ierr != 0 || gpar[2] <= 0.0 )
1637  {
1638  found = NOCONV;
1639  istat = 3;
1640  }
1641  else if ( (oldchi - chisqr)/ chisqr < GCHIFND )
1642  found = TRUE;
1643  }
1644 
1645  /* Is the source still in the image and has it a resonable shape? */
1646 
1647  lamda = 0.0;
1648  ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma);
1649  if ( ierr != 0 )
1650  {
1651  found = NOCONV;
1652  istat = 3;
1653  }
1654  else
1655  {
1656  indx = CGN_NINT(xypos[0]);
1657  sumi = (float) (gpar[1] + image[2]);
1658  indy = CGN_NINT(sumi);
1659  if ( indy < 0 || indy >= npix[1] )
1660  {
1661  found = OUTSIDE;
1662  istat = 2;
1663  }
1664  else
1665  indx += (*npix) * indy;
1666  }
1667  (void) cpl_free( (char *) xpos );
1668  (void) cpl_free( (char *) yfit );
1669  (void) cpl_free( (char *) gfit );
1670 
1671  if ( found == TRUE )
1672  {
1673  xypos[1] = sumi;
1674  xysig[1] = (float) gpar[2];
1675  xyerr[1] = (float) sqrt( sigma * chisqr );
1676  *xyval = p_img[indx];
1677  }
1678  }
1679  (void) cpl_free( ( char *) kry );
1680  }
1681  }
1682 
1683 return istat;
1684 }
static int CGN_NINT(float a)
finds absolute value of nearest integer
static void Crhox(float *p_img, int *npix, int *image, int *lnew, double *krx)
compute X-marginal vector KRX.
static int LSQFIT(double *xdat, double *data, int ndim, double *gpar, float *lamda, double *dfit, double *chisqr, double *sigma)
least squares fit to a non-linear function
int uves_physmod_stacen(float *p_img, int dimx, int dimy, char meth, int *image, float *xout, float *yout, float *xerr, float *yerr, float *xsig, float *ysig, float *xyval, int *stat)
Routines used to do Gaussian fit to a line.
static float FCHIS(double *data, int ndim, int nfree, int mode, double *dfit)
evaluate reduced chi square for fit to data
static double ERFCC(double xx)
returns complementary error function EFC( xx )
int uves_physmod_cstacen(char meth, float *p_img, int *npix, int *image, float *xypos, float *xyerr, float *xysig, float *xyval)
Routines used to do Gaussian fit to a line.
static double GAUSFU(double xx, double *gpar)
static int Cserch(double *marg, int ndim, int ign, int *lmin, int *lmax, float *s_cent, float *s_width)
search a star from a marginal distribution
static int MATINV(double(*matrix)[MAXPAR], int nfree)
calculate the inverse of a matrix
static void GAUSDE(double xdat, double *gpar, double *deriv)
evaluates derivatives of function for least squares search with shape of a gaussian distribution ...
static void Crhoy(float *p_img, int *npix, int *image, int *lnew, double *kry)
compute Y-marginal vector KRY.
static int Ckapsig(float *val, int nval, int iter, float akap, float *cons, float *rms, int *npts)
selects a constant mean through a kap*sig clipping