UVES Pipeline Reference Manual  5.4.0
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  * $Log: not supported by cvs2svn $
26  * Revision 1.9 2007/06/06 08:17:33 amodigli
27  * replace tab with 4 spaces
28  *
29  * Revision 1.8 2006/08/23 15:41:06 amodigli
30  * removed warning from checks on line length
31  *
32  * Revision 1.7 2006/08/11 14:56:05 amodigli
33  * removed Doxygen warnings
34  *
35  * Revision 1.6 2006/06/20 10:56:56 amodigli
36  * cleaned output, added units
37  *
38  * Revision 1.5 2006/06/20 08:25:56 amodigli
39  * fixed doxigen warnings
40  *
41  * Revision 1.4 2006/06/13 11:59:51 jmlarsen
42  * Fixed doc. bug
43  *
44  * Revision 1.3 2006/06/08 11:01:50 amodigli
45  * fixed some warnings
46  *
47  * Revision 1.2 2006/06/01 14:43:17 jmlarsen
48  * Added missing documentation
49  *
50  * Revision 1.1 2006/02/03 07:46:30 jmlarsen
51  * Moved recipe implementations to ./uves directory
52  *
53  * Revision 1.13 2006/01/20 10:05:49 jmlarsen
54  * Inserted missing doxygen end tag
55  *
56  * Revision 1.12 2006/01/13 09:54:42 amodigli
57  * Fixed some bugs: improved agreement with MIDAS version
58  *
59  * Revision 1.11 2006/01/09 14:05:42 amodigli
60  * Fixed doxigen warnings
61  *
62  * Revision 1.10 2006/01/05 14:29:59 jmlarsen
63  * Removed newline characters from output strings
64  *
65  * Revision 1.9 2005/12/20 08:11:44 jmlarsen
66  * Added CVS entry
67  *
68  */
69 
70 #ifdef HAVE_CONFIG_H
71 # include <config.h>
72 #endif
73 
74 /*---------------------------------------------------------------------------*/
81 /*---------------------------------------------------------------------------*/
82 
83 /* code derived by MIDAS cstacen.c */
84 /*-----------------------------------------------------------------------------
85  Includes
86  ----------------------------------------------------------------------------*/
87 /* definition of the used functions in this module */
88 
89 #include "uves_physmod_cstacen.h"
90 #include <cpl.h>
91 #include <math.h>
92 /*-----------------------------------------------------------------------------
93  Defines
94  ----------------------------------------------------------------------------*/
95 /* Define _POSIX_SOURCE to indicate that this is a POSIX program */
96 /* replaced osmmget by cpl_calloc */
97 /* replaced osmmfree by cpl_free */
98 #define _POSIX_SOURCE 1
99 
100 /* define some macros and constants */
101 
102 #ifndef PI
103 #define PI 3.14159265358979325e0
104 #endif
105 /*
106 #ifndef true
107 #define true 1
108 #define false 0
109 #endif
110 */
111 
112 /* Constants used by the moment centering */
113 
114 #define MINVAL 1.0e-37
115 #define MMXITER 8
116 #define SMALL 1.0e-20
117 
118 
119 /* Constants used by the gaussian centering */
120 
121 #define MAXPAR 4
122 #define IGNORE 2
123 #define NOCONV -1
124 #define OUTSIDE -2
125 #define GMXITER 50
126 #define GCHIMAX 5.0e+16
127 #define GCHIFND 0.005
128 
129 
130 #define MYMIN(a,b) ((a) > (b) ? (b) : (a))
131 #define MYMAX(a,b) ((b) > (a) ? (b) : (a))
132 
135 /*-----------------------------------------------------------------------------
136  Functions prototypes
137  ----------------------------------------------------------------------------*/
138 
139 /*-----------------------------------------------------------------------------
140  Static variables
141  ----------------------------------------------------------------------------*/
142 /*-----------------------------------------------------------------------------
143  Functions code
144  ----------------------------------------------------------------------------*/
145 
146 /*---------------------------------------------------------------------------*/
185 /*---------------------------------------------------------------------------*/
186 
187 int
188 uves_physmod_stacen(float* p_img, int dimx, int dimy, char meth, int* image,
189  float* xout, float* yout, float* xerr, float* yerr,
190  float* xsig, float* ysig, float* xyval, int* stat )
191 
192 {
193 int npix[2], imap[4];
194 float xypos[2], xyerr[2], xysig[2];
195 
196 npix[0] = dimx;
197 npix[1] = dimy;
198 
199 imap[0] = image[0] - 1;
200 imap[1] = image[1] - 1;
201 imap[2] = image[2] - 1;
202 imap[3] = image[3] - 1;
203 
204 /*
205  uves_msg("Input=: npix[0]=%d npix[1]=%d",npix[0],npix[1]);
206  uves_msg("Input=: imap[0]=%d imap[1]=%d imap[2]=%d imap[3]=%d",imap[0],imap[1],imap[2],imap[3]);
207 */
208 
209 *stat = uves_physmod_cstacen(meth, p_img, npix, imap, xypos, xyerr, xysig, xyval );
210 
211 *xout = xypos[0] + 1;
212 *yout = xypos[1] + 1;
213 
214 
215 *xerr = xyerr[0];
216 *yerr = xyerr[1];
217 *xsig = xysig[0];
218 *ysig = xysig[1];
219 /*
220  uves_msg("xout=%f,yout=%f,xerr=%f,yerr=%f,xsig=%f,ysig=%f",
221  *xout,*yout,*xerr,*yerr,*xsig,*ysig);
222  */
223 
224  return 0;
225 }
226 /*---------------------------------------------------------------------------*/
232 static int CGN_NINT(float a){
233  int res=0;
234 
235  res = (a) > 0 ? floor(a+0.5) : -floor(a-0.5);
236  return res;
237 }
238 
239 /*---------------------------------------------------------------------------*/
266 /*---------------------------------------------------------------------------*/
267 
268 
269 static int Ckapsig( float *val, int nval, int iter, float akap,
270  float *cons, float *rms, int *npts )
271 
272 {
273 
274 register int ii=0;
275 register int it=0;
276 register int nr=0;
277 
278 int nr_old=0;
279 
280 float clip=0;
281 float dels=0;
282 float delv=0;
283 float mean=0;
284 float msq=0;
285 float sum=0;
286 float *vsq=NULL;
287 
288 if ( nval < 2 ) return (-1);
289 
290 
291 /* initialize mean value */
292 
293 mean = 0.0;
294 for (ii=0; ii<nval; ii++) mean += val[ii];
295 mean /= (float) (nval);
296 msq = mean * mean;
297 
298 
299 /* initialize RMS */
300 
301 vsq = (float *) cpl_calloc( nval, sizeof( float ));
302 
303 dels = 0.0;
304 for (ii=0; ii<nval; ii++)
305  {
306  vsq[ii] = val[ii] * val[ii];
307  delv = MYMAX( 0.0, vsq[ii] + msq - (2.0 * mean * val[ii]));
308  dels += delv;
309  }
310 
311 *rms = (float) sqrt( MYMAX( MINVAL, dels / (nval-1)));
312 clip = akap * (*rms);
313 
314 
315 /* iterate */
316 
317 nr_old = 0;
318 for (it=0; it<iter; it++)
319  {
320  nr = 0;
321  sum = 0.0;
322  dels = 0.0;
323 
324  for ( ii = 0; ii < nval; ii++ )
325  {
326  if ( fabs( val[ii] - mean ) < (double) clip )
327  {
328  nr++;
329  delv = MYMAX( 0.0, vsq[ii] + msq - 2.0 * mean * val[ii]);
330  dels += delv;
331  sum += val[ii];
332  }
333  }
334 
335  if ( nr <= 2 || nr == nr_old ) goto end_of_it;
336 
337 
338  /* define new rms and mean value */
339 
340  nr_old = nr;
341  *rms = (float) sqrt( MYMAX( MINVAL, dels / (nr-1)));
342  clip = akap * *rms;
343  mean = sum / nr_old;
344  msq = mean * mean;
345  }
346 
347 end_of_it:
348 *cons = mean; /* exit */
349 *npts = nr;
350 
351 cpl_free(vsq);
352 return (0);
353 }
354 
355 
356 /*---------------------------------------------------------------------------*/
371 /*---------------------------------------------------------------------------*/
372 
373 static int MATINV( double (*matrix)[MAXPAR], int nfree )
374 {
375 
376 register int ii=0;
377 register int jj=0;
378 register int kk=0;
379 
380 int evin=0;
381 int row=0;
382 int per[MAXPAR];
383 double even=0.;
384 double mjk=0.;
385 double rowmax=0.;
386 double hv[MAXPAR];
387 
388 
389 for ( ii = 0; ii < nfree; ii++ ) per[ii] = ii; /* set permutation array */
390 
391 for ( jj = 0; jj < nfree; jj++ ) /* in j-th column, ... */
392  {
393  rowmax = fabs( matrix[jj][jj] ); /* determine row with ... */
394  row = jj; /* largest element. */
395  for ( ii = jj + 1; ii < nfree; ii++ )
396  {
397  if ( fabs( matrix[ii][jj] ) > rowmax )
398  {
399  rowmax = fabs( matrix[ii][jj] );
400  row = ii;
401  }
402  }
403 
404  if (fabs(matrix[row][jj]) < SMALL) /* determinant is zero! */
405  return (1);
406 
407  if ( row > jj ) /* if largest element not ...*/
408  {
409  for ( kk = 0; kk < nfree; kk++ ) /* on diagonal, then ... */
410  {
411  even = matrix[jj][kk]; /* permutate rows. */
412  matrix[jj][kk] = matrix[row][kk];
413  matrix[row][kk] = even;
414  }
415  evin = per[jj]; /* keep track of permutation */
416  per[jj] = per[row];
417  per[row] = evin;
418  }
419 
420  even = 1.0 / matrix[jj][jj]; /* modify column */
421  for (ii=0; ii<nfree; ii++)
422  matrix[ii][jj] *= even;
423  matrix[jj][jj] = even;
424  for (kk=0; kk<jj; kk++)
425  {
426  mjk = matrix[jj][kk];
427  for ( ii = 0; ii < jj; ii++ )
428  matrix[ii][kk] -= matrix[ii][jj] * mjk;
429  for ( ii = jj + 1; ii < nfree; ii++ )
430  matrix[ii][kk] -= matrix[ii][jj] * mjk;
431  matrix[jj][kk] = -even * mjk;
432  }
433 
434  for ( kk = jj + 1; kk < nfree; kk++ )
435  {
436  mjk = matrix[jj][kk];
437  for ( ii = 0; ii < jj; ii++ )
438  matrix[ii][kk] -= matrix[ii][jj] * mjk;
439  for ( ii = jj + 1; ii < nfree; ii++ )
440  matrix[ii][kk] -= matrix[ii][jj] * mjk;
441  matrix[jj][kk] = -even * mjk;
442  }
443  }
444 
445 for ( ii = 0; ii < nfree; ii++ ) /* finally, repermute the ...*/
446  {
447  for ( kk = 0; kk < nfree; kk++ ) /* columns. */
448  {
449  hv[per[kk]] = matrix[ii][kk];
450  }
451  for ( kk = 0; kk < nfree; kk++)
452  {
453  matrix[ii][kk] = hv[kk];
454  }
455  }
456 
457 return 0;
458 }
459 
460 /*---------------------------------------------------------------------------*/
471 /*---------------------------------------------------------------------------*/
472 
473 static double ERFCC( double xx )
474 
475 {
476 
477 double t=0.;
478 double z=0.;
479 double ans=0.;
480 double zz=0.;
481 double zozo=0.;
482 
483 
484 z = fabs( xx );
485 t = 1.0 / (1.0 + (0.5 * z));
486 
487 /*
488 ans = t * exp( -z * z - 1.26551223 + t * ( 1.00002368 + t *
489  ( 0.37409196 + t * ( 0.09678418 + t *
490  (-0.18628806 + t * ( 0.27886807 + t *
491  (-1.13520398 + t * ( 1.48851587 + t *
492  (-0.82215223 + t * 0.17087277 )))))))));
493 
494 
495 the original code above didn't work on Red Hat Linux 5.2
496 with CENTER/GAUSS where the main program is Fortran code
497 neither on Alpha nor on Intel PC (using f2c) ...
498 
499 however it works on SUSE Linux 6.xx on Intel (using g77)
500 
501 therefore this work around which seemed to solve the problem, KB 000522 */
502 
503 zz = -z * z - 1.26551223 + t * ( 1.00002368 + t *
504  ( 0.37409196 + t * ( 0.09678418 + t *
505  (-0.18628806 + t * ( 0.27886807 + t *
506  (-1.13520398 + t * ( 1.48851587 + t *
507  (-0.82215223 + t * 0.17087277 ))))))));
508 
509 
510 if (zz < -500.0)
511  zozo = 0.0;
512 else
513  zozo = exp(zz);
514 
515 ans = t * zozo;
516 
517 
518 return (xx >= 0.0 ? ans : 2.0 - ans);
519 }
520 
521 /*---------------------------------------------------------------------------*/
531 /*---------------------------------------------------------------------------*/
532 
533 static double GAUSFU( double xx, double *gpar )
534 {
535 
536 double rc=0.;
537 double dd=0.;
538 
539 static int init = TRUE;
540 static double sqrt_2;
541 static double rc1;
542 
543 
544 if ( init )
545  {
546  sqrt_2 = sqrt( 2.0 );
547  rc1 = sqrt(PI)/sqrt_2;
548  init = FALSE;
549  }
550 
551 rc = 1.0 / (sqrt_2 * gpar[2]);
552 dd = xx - gpar[1];
553 dd = ERFCC(rc * (dd - 0.5)) - ERFCC(rc * (dd + 0.5));
554 return ( gpar[3] + rc1 * gpar[0] * gpar[2] * dd );
555 }
556 
557 /*---------------------------------------------------------------------------*/
569 /*---------------------------------------------------------------------------*/
570 
571 static void GAUSDE( double xdat, double *gpar, double *deriv )
572 {
573 double temp=0.;
574 double tempp=0.;
575 double x1=0.;
576 double x2=0.;
577 double zz=0.;
578 double zx=0.;
579 double dv1=0.;
580 double dv2=0.;
581 static double sqrt_2;
582 
583 register int jj=0;
584 static int init = TRUE;
585 
586 
587 if ( init )
588  {
589  sqrt_2 = sqrt( 2.0 );
590  init = FALSE;
591  }
592 
593 temp = sqrt_2 * gpar[2];
594 tempp = xdat - gpar[1];
595 x1 = (tempp - 0.5) / temp;
596 x2 = (tempp + 0.5) / temp;
597 zz = tempp / gpar[2] ;
598 
599 if ( ((zz * zz) - 50.0) < 0.0 )
600  {
601  deriv[0] = (GAUSFU( xdat, gpar ) - gpar[3]) / gpar[0];
602 
603  zx = (-x1) * x1;
604  if ( zx < -200.0) /* zx always < 0 */
605  dv1 = 0.0; /* e**(-200) is = 0.0 ... */
606  else
607  dv1 = exp(zx);
608  zx = (-x2) * x2;
609  if ( zx < -200.0)
610  dv2 = dv1;
611  else
612  dv2 = dv1 - exp(zx);
613  deriv[1] = gpar[0] * dv2;
614 
615  /* for (x1 * x1) > 400 we got floating point exceptions on DEC Alpha
616  deriv[1] = gpar[0] * (exp( -x1 * x1 ) - exp( -x2 * x2 ));
617  */
618 
619  deriv[2] = deriv[1] * zz;
620  }
621 else
622  for (jj=0; jj<3; jj++) deriv[jj] = 0.0;
623 
624 deriv[3] = 1.0;
625 }
626 
627 /*---------------------------------------------------------------------------*/
643 /*---------------------------------------------------------------------------*/
644 
645 static float FCHIS(double *data,int ndim,int nfree,int mode,double *dfit) {
646 
647 register int ii=0;
648 
649 double diff=0.;
650 double weight=0.;
651 double chisq=0.;
652 
653 
654 if ( nfree > 0 )
655  {
656  chisq = 0.0;
657 
658  for (ii=0; ii<ndim; ii++)
659  {
660  if ( mode < 0 )
661  {
662  if ( *data < 0 )
663  weight = -1. / *data;
664  else if ( *data == 0 )
665  weight = 1.0;
666  else
667  weight = 1. / *data;
668  }
669  else
670  weight = 1.0;
671 
672  diff = (*data) - (*dfit);
673  data++; dfit++;
674  chisq += weight * diff * diff;
675  }
676  return (chisq / nfree);
677  }
678 
679 else
680  return 0.0;
681 }
682 
683 /*---------------------------------------------------------------------------*/
711 static int LSQFIT( double *xdat, double *data, int ndim,
712  double *gpar, float *lamda, double *dfit,
713  double *chisqr, double *sigma )
714 
715 {
716 
717 register int icnt=0;
718 register int ii=0;
719 register int jj=0;
720 register int kk=0;
721 
722 int nfree=0;
723 
724 double chisq1=0;
725 double b[MAXPAR], beta[MAXPAR], deriv[MAXPAR],
726  array[MAXPAR][MAXPAR], alpha[MAXPAR][MAXPAR];
727 
728 
729 nfree = ndim - MAXPAR;
730 *sigma = 0.0;
731 if ( nfree < 1 || fabs( (double) *gpar ) < SMALL ) return (1);
732 
733 
734 /* evaluate ALPHA and BETA matrices */
735 
736 for (ii=0; ii<MAXPAR; ii++)
737  {
738  beta[ii] = 0.0;
739  for (jj=0; jj<=ii; jj++) alpha[ii][jj] = 0.0;
740  }
741 
742 for (ii=0; ii<ndim; ii++)
743  {
744  GAUSDE( xdat[ii], gpar, deriv ); /* here we divide by gpar[1] */
745 
746  for (jj=0; jj<MAXPAR; jj++)
747  {
748  beta[jj] += (data[ii] - GAUSFU( xdat[ii], gpar )) * deriv[jj];
749  for (kk=0; kk<=jj; kk++)
750  alpha[jj][kk] += deriv[jj] * deriv[kk];
751  }
752  }
753 
754 for (ii=0; ii<MAXPAR; ii++)
755  {
756  for (jj=0; jj<=ii; jj++)
757  alpha[jj][ii] = alpha[ii][jj];
758  }
759 
760 
761 /* invert matrix */
762 
763 if ( *lamda < SMALL)
764  {
765  if (MATINV(alpha,MAXPAR) == 1) return (2); /* determinant -> 0.0 */
766 
767  *sigma = MYMAX( 0.0, alpha[1][1] );
768  }
769 
770 else /* evaluate chi square at starting point */
771  {
772  for (ii=0; ii<ndim; ii++)
773  dfit[ii] = GAUSFU( xdat[ii], gpar );
774 
775  chisq1 = FCHIS( data, ndim, nfree, 0, dfit );
776 
777  icnt = 0; /* invert matrix */
778 loop:
779  for ( jj = 0; jj < MAXPAR; jj++ )
780  {
781  for ( kk = 0; kk < MAXPAR; kk++ )
782  {
783  if (fabs( alpha[jj][jj] ) < 1.e-15 || fabs( alpha[kk][kk] ) < 1.e-15)
784  return 2;
785  array[jj][kk] = alpha[jj][kk] /
786  sqrt( alpha[jj][jj] * alpha[kk][kk] ) ;
787  }
788  array[jj][jj] = 1.0 + *lamda;
789  }
790 
791  (void) MATINV( array, MAXPAR );
792 
793  for ( jj = 0; jj < MAXPAR; jj++ )
794  {
795  b[jj] = gpar[jj] ;
796  for ( kk = 0; kk < MAXPAR ; kk++ )
797  {
798  b[jj] += beta[kk] * array[jj][kk] /
799  sqrt( alpha[jj][jj] * alpha[kk][kk] );
800  }
801  }
802 
803 /* if chi square increased, increase LAMDA and try again */
804 
805  for (ii=0; ii<ndim; ii++)
806  dfit[ii] = GAUSFU( xdat[ii], b );
807 
808  *chisqr = FCHIS( data, ndim, nfree, 0, dfit );
809 
810  if ( chisq1 - *chisqr < 0.0 )
811  {
812  if (++icnt < 60)
813  {
814  *lamda *= 10;
815  goto loop;
816  }
817  else
818  return (2);
819  }
820 
821  for (jj=0; jj<MAXPAR; jj++) gpar[jj] = b[jj];
822  *lamda /= 10.0;
823  }
824 
825 return 0;
826 }
827 
828 
829 /*---------------------------------------------------------------------------*/
849 static void Crhox( float *p_img, int *npix, int *image,
850  int *lnew, double *krx )
851 
852 {
853 register int nxdim=0;
854 register int ix=0;
855 register int iy=0;
856 
857 int nrx=0;
858 int nry=0;
859 
860 double sum=0;
861 
862 
863 
864 
865 /* original FORTRAN code:
866 
867  IXA = IMAP(1) IMAP(1-4) => image[0-3]
868  IXE = IMAP(2)
869  IYA = IMAP(3)
870  JYA = IYA + JY - 1 JY => lnew[0]
871  JYE = IYA + LY - 1 LY => lnew[1]
872  M = 1
873 
874  DO 200 J=IXA,IXE
875  ISUM = 0.D0
876  DO 100 K=JYA,JYE
877  ISUM = ISUM + AIMG(J,K)
878 100 CONTINUE
879  KRX(M) = ISUM
880  M = M + 1
881 200 CONTINUE
882 */
883 
884 
885 nrx = image[1] - image[0] + 1;
886 nry = lnew[1] - lnew[0] + 1;
887 nxdim = *npix;
888 p_img += nxdim * (image[2] + lnew[0]);
889 
890 for (ix=0; ix<nrx; ix++)
891  {
892  sum = 0.0;
893  for (iy=0; iy<nry*nxdim; iy+=nxdim) sum += p_img[iy];
894  p_img++;
895  *krx++ = sum;
896  }
897 }
898 
899 
900 
901 /*---------------------------------------------------------------------------*/
921 static void Crhoy( float *p_img, int *npix, int *image,
922  int *lnew, double *kry )
923 
924 {
925 register int nxdim=0;
926 register int ix=0;
927 register int iy=0;
928 
929 int nrx=0;
930 int nry=0;
931 
932 double sum=0;
933 
934 
935 
936 /* original FORTRAN code:
937 
938  IXA = IMAP(1) IMAP(1-4) => image[0-3]
939  IYA = IMAP(3)
940  IYE = IMAP(4)
941  JXA = IXA + JX - 1 JX => lnew[0]
942  JXE = IXA + LX - 1 LX => lnew[1]
943 
944  M = 1
945  DO 200 J=IYA,IYE
946  ISUM = 0.D0
947  DO 100 K=JXA,JXE
948  ISUM = ISUM + AIMG(K,J)
949 100 CONTINUE
950  KRY(M) = ISUM
951  M = M + 1
952 200 CONTINUE
953 
954 */
955 
956 nrx = lnew[1] - lnew[0] + 1;
957 nry = image[3] - image[2] + 1;
958 nxdim = *npix;
959 p_img += (nxdim * image[2]) + (image[0] + lnew[0]);
960 
961 for (iy=0; iy<nry; iy++)
962  {
963  sum = 0.0;
964  for (ix=0; ix<nrx; ix++) sum += p_img[ix];
965  p_img += nxdim;
966  *kry++ = sum;
967  }
968 }
969 
970 /*---------------------------------------------------------------------------*/
994 static int Cserch( double *marg, int ndim, int ign,
995  int *lmin, int *lmax, float *s_cent, float *s_width )
996 
997 {
998 
999 register int ii=0;
1000 
1001 int ql=0;
1002 int ibgn=0;
1003 int icrowd=0;
1004 int iend=0;
1005 int imax=0;
1006 int imin=0;
1007 int indx=0;
1008 
1009 double dxk=0.;
1010 double diff=0.;
1011 double drmn=0.;
1012 double drmx=0.;
1013 double sum=0.;
1014 double *work=NULL;
1015 
1016 
1017 ibgn = ign;
1018 iend = ndim - ign -1; /* ojo */
1019 
1020 /* create workspace to store the derivative of the marginal data */
1021 
1022 work = (double *) cpl_calloc( ndim , sizeof( double ));
1023 
1024 
1025 /* find maximum and minimum derivative of MARG */
1026 
1027 imin = imax = 0;
1028 drmn = drmx = 0.0;
1029 for (ii = ibgn; ii < iend; ii++ )
1030  {
1031  diff = marg[ii+1] - marg[ii-1];
1032  work[ii] = marg[ii+2] - marg[ii-2] + (2 * diff);
1033  if ( work[ii] >= drmx )
1034  {
1035  drmx = work[ii];
1036  imax = ii;
1037  }
1038  if (
1039  work[ii] <= drmn )
1040  { drmn = work[ii];
1041  imin = ii;
1042  }
1043  }
1044 
1045 
1046 /* crowded ? */
1047 
1048 icrowd = 0;
1049 if (imin <= imax ) /* bright source to the left, compute right a new minima */
1050  {
1051  if ( ndim - imax > imin )
1052  {
1053  icrowd = -1;
1054  drmn = drmx;
1055  for ( ii = imax+1; ii < iend; ii++ )
1056  {
1057  if ( work[ii] < drmn )
1058  {
1059  drmn = work[ii];
1060  imin = ii;
1061  }
1062  }
1063  }
1064  else /* bright source to the right, compute left a new maxima */
1065  {
1066  icrowd = 1;
1067  drmx = drmn;
1068  for ( ii = ibgn; ii < imin; ii++ )
1069  {
1070  if ( work[ii] >= drmx )
1071  {
1072  drmx = work[ii];
1073  imax = ii;
1074  }
1075  }
1076  }
1077  }
1078 
1079 
1080 /* compute estimates of image centre and width */
1081 
1082 *s_cent = ((float)(imax + imin)) * 0.5;
1083 *s_width = imin - imax;
1084 
1085 sum = 0.0;
1086 for ( ii = imax; ii <= imin; ii++ ) sum += work[ii];
1087 
1088 diff = drmx - drmn;
1089 if ( fabs(diff) > SMALL)
1090  {
1091  dxk = sum * *s_width / ( (*s_width+1.0)*diff );
1092  *s_cent += dxk;
1093  }
1094 *s_width /= 2;
1095 indx = CGN_NINT(*s_cent);
1096 if (indx < 0)
1097  {
1098  *s_cent = 0.0;
1099  indx = 0;
1100  }
1101 else if (indx >= ndim)
1102  {
1103  *s_cent = (float)(ndim-1);
1104  indx = ndim-1;
1105  }
1106 
1107 
1108 /* find low- (left-) side local minimum */
1109 
1110 ql = indx - 2;
1111 *lmin = 0;
1112 if (ql < 2) goto next_step;
1113 
1114 low_loop:
1115 ql --;
1116 if (ql <= 0 ) goto next_step;
1117 if (ql == 1) goto lo5;
1118 if (ql == 2) goto lo4;
1119 if (ql == 3) goto lo3;
1120 
1121 if (marg[ql] > marg[ql-4]) goto low_loop;
1122 lo3:
1123 if (marg[ql] > marg[ql-3]) goto low_loop;
1124 lo4:
1125 if (marg[ql] > marg[ql-2]) goto low_loop;
1126 lo5:
1127 if (marg[ql] > marg[ql-1]) goto low_loop;
1128 
1129 *lmin = ql + 1;
1130 
1131 
1132 /* find high- (right-) side local minimum */
1133 
1134 next_step:
1135 ql = indx + 2;
1136 *lmax = ndim - 1;
1137 ii = ndim - ql;
1138 if (ii < 3) goto end_of_it;
1139 
1140 hi_loop:
1141 ql ++ ;
1142 ii = ndim - ql;
1143 if (ii == 1 ) goto end_of_it;
1144 if (ii == 2) goto hi5;
1145 if (ii == 3) goto hi4;
1146 if (ii == 4) goto hi3;
1147 
1148 if (marg[ql] > marg[ql+4]) goto hi_loop;
1149 hi3:
1150 if (marg[ql] > marg[ql+3]) goto hi_loop;
1151 hi4:
1152 if (marg[ql] > marg[ql+2]) goto hi_loop;
1153 hi5:
1154 if (marg[ql] > marg[ql+1]) goto hi_loop;
1155 *lmax = ql - 1;
1156 
1157 end_of_it:
1158 (void) cpl_free( (char *) work );
1159 return (icrowd);
1160 
1161 }
1162 
1163 /*---------------------------------------------------------------------------*/
1201 /*---------------------------------------------------------------------------*/
1202 
1203 /* ------------------------------------------*/
1204 /* here starts the code of the main function */
1205 /* ------------------------------------------*/
1206 
1207 
1208 int uves_physmod_cstacen(char meth, float* p_img, int* npix, int* image,
1209  float* xypos, float* xyerr, float* xysig, float* xyval )
1210 
1211 {
1212 
1213 register int it=0;
1214 register int ix=0;
1215 register int iy=0;
1216 
1217 int bgnr=0;
1218 int indx=0;
1219 int indy=0;
1220 int istat=0;
1221 int nrx=0;
1222 int nry=0;
1223 int nval=0;
1224 int ifram[4];
1225 
1226 float bgval=0.;
1227 float clip=0.;
1228 float rms=0.;
1229 float xmom=0.;
1230 float ymom=0.;
1231 float source=0.;
1232 float sumi=0.;
1233 float xold=0.;
1234 float yold=0.;
1235 
1236 float *p_bgn=NULL;
1237 float *p_edge=NULL;
1238 
1239 
1240 istat = 0; /* initialize */
1241 p_bgn = p_img;
1242 for (ix=0; ix<4; ix++)
1243  ifram[ix] = image[ix] + 1; /* 1 ---> ndim */
1244 nrx = ifram[1] - ifram[0] + 1;
1245 nry = ifram[3] - ifram[2] + 1;
1246 
1247 xypos[0] = (ifram[0] + ifram[1]) * 0.5; /* init to center pixel */
1248 xypos[1] = (ifram[2] + ifram[3]) * 0.5;
1249 xyerr[0] = xyerr[1] = 0.0;
1250 xysig[0] = xysig[1] = 0.0;
1251 *xyval = 0.0;
1252 
1253 /* MOMENT centering */
1254 
1255 if ( meth != 'G' && meth != 'g' )
1256  {
1257  int kk, istr, iend;
1258 
1259  xold = yold = -1.0; /* find bgval and rms from edge pixels */
1260 
1261  p_img += (ifram[0] - 1) + (npix[0] * (ifram[2] - 1));
1262 
1263  /* collect edge pixels */
1264 
1265  if (nry > 1)
1266  {
1267  nval = (2 * nrx) + (2 * (nry-2));
1268  p_edge = (float *) cpl_calloc( nval , sizeof( float ));
1269 
1270  for (ix=0; ix<nrx;ix++)
1271  *p_edge++ = p_img[ix];
1272 
1273  p_img += *npix;
1274  for (iy=0; iy<(nry-2); iy++)
1275  {
1276  *p_edge++ = p_img[0];
1277  *p_edge++ = p_img[nrx - 1];
1278  p_img += npix[0];
1279  }
1280  for (ix=0; ix<nrx; ix++) *p_edge++ = p_img[ix];
1281  }
1282  else
1283  {
1284  nval = nrx;
1285  p_edge = (float *) cpl_calloc( nval , sizeof( float ));
1286 
1287  for (ix=0; ix<nrx;ix++)
1288  *p_edge++ = p_img[ix];
1289  }
1290 
1291  p_img = p_bgn;
1292 
1293  p_edge -= nval;
1294  (void) Ckapsig( p_edge, nval, 5, 2.0, &bgval, &rms, &bgnr );
1295  (void) cpl_free( (char *) p_edge );
1296 
1297 
1298  /* calculate moment for pixel values > 3 * RMS above BGVAL */
1299 
1300  clip = 3.0 * rms;
1301 
1302  for (it=0; it<MMXITER; it++) /* iteration loop */
1303  {
1304  sumi = xmom = ymom = 0.0;
1305  p_img += ifram[0] - 1 + (npix[0] * (ifram[2] - 1));
1306  for (iy=0; iy<nry; iy++)
1307  {
1308  for (ix=0; ix<nrx; ix++)
1309  {
1310  if ( (source = p_img[ix] - bgval) > clip )
1311  {
1312  sumi += source;
1313  xmom += source * (ifram[0] + ix);
1314  ymom += source * (ifram[2] + iy);
1315  }
1316  }
1317  p_img += npix[0];
1318  }
1319  p_img = p_bgn; /* reset to start of array */
1320 
1321  if ((nrx < 3) || (nry < 3))
1322  {
1323  xysig[0] = nrx;
1324  xysig[1] = nry;
1325  istat = 1;
1326  if ( sumi > 0.0 )
1327  {
1328  xypos[0] = xmom / sumi;
1329  xypos[1] = ymom / sumi;
1330  }
1331  else
1332  {
1333  istat = 2;
1334  xypos[0] = (ifram[0] + ifram[1]) * 0.5;
1335  xypos[1] = (ifram[2] + ifram[3]) * 0.5;
1336  }
1337  indx = CGN_NINT(xypos[0]-1);
1338  indy = CGN_NINT(xypos[1]-1);
1339  *xyval = p_img[indx + ((*npix) * indy)];
1340  goto end_of_iter; /* EXIT iteration loop */
1341  }
1342 
1343  if (sumi > 0.0) /* only positive sources */
1344  {
1345  xypos[0] = xmom / sumi;
1346  xypos[1] = ymom / sumi;
1347  xysig[0] = nrx;
1348  xysig[1] = nry;
1349 
1350  if ( xold == xypos[0] && yold == xypos[1] )
1351  {
1352  int nr = 0;
1353  double xdif, ydif, xrms, yrms;
1354 
1355  xrms = yrms = sumi = 0.0;
1356  p_img += ifram[0] - 1 + (npix[0] * (ifram[2] - 1));
1357  for (iy=0; iy<nry; iy++ )
1358  {
1359  for (ix=0; ix<nrx; ix++)
1360  {
1361  if ( (source = p_img[ix] - bgval) > clip )
1362  {
1363  xdif = (ifram[0] + ix) - xypos[0];
1364  ydif = (ifram[2] + iy) - xypos[1];
1365  xrms += fabs( source * xdif *xdif );
1366  yrms += fabs( source * ydif *ydif );
1367  sumi += source;
1368  nr++;
1369  }
1370  }
1371  p_img += npix[0];
1372  }
1373  p_img = p_bgn;
1374 
1375  indx = CGN_NINT(xypos[0]-1) + (npix[0] * CGN_NINT(xypos[1]-1));
1376  *xyval = p_img[indx];
1377  xysig[0] = (float) sqrt(xrms /(sumi+ *xyval - bgval));
1378  xysig[1] = (float) sqrt(yrms /(sumi+ *xyval - bgval));
1379  xyerr[0] = (float) (xysig[0] / sqrt( (double) (nr - 1)));
1380  xyerr[1] = (float) (xysig[1] / sqrt( (double) (nr - 1)));
1381  goto end_of_iter; /* succesful return */
1382  }
1383 
1384 
1385  xold = xypos[0];
1386  yold = xypos[1];
1387  }
1388  else
1389  {
1390  istat = 2;
1391  xypos[0] = (ifram[0] + ifram[1]) * 0.5;
1392  xypos[1] = (ifram[2] + ifram[3]) * 0.5;
1393  indx = CGN_NINT(xypos[0]-1);
1394  indy = CGN_NINT(xypos[1]-1);
1395  *xyval = p_img[indx + ((*npix) * indy)];
1396  goto end_of_iter; /* EXIT iteration loop */
1397  }
1398 
1399 
1400  /* crowded or weak source conditions */
1401 
1402  indx = CGN_NINT(xypos[0]-1) + (npix[0] * CGN_NINT(xypos[1]-1));
1403  if ( (*xyval = p_img[indx] - bgval) <= clip )
1404  {
1405  xysig[0] = xysig[1] = 0.0;
1406  istat = 1;
1407  goto end_of_iter; /* EXIT iteration loop */
1408  }
1409 
1410 
1411  /* find extent of source i.e. delete spikes, etc. */
1412 
1413  ix = CGN_NINT( xypos[0] ); /* ix, iy = 1,2,... */
1414  iy = CGN_NINT( xypos[1] );
1415  kk = npix[0] * (iy - 1);
1416  istr = ifram[0];
1417  source = p_img[ix-1 + kk] - bgval;
1418  while ( source > clip && ix >= istr )
1419  {
1420  ifram[0] = ix;
1421  source = p_img[ix-1 + kk] - bgval;
1422  ix --;
1423  }
1424 
1425  ix = CGN_NINT( xypos[0] );
1426  iend = ifram[1];
1427  source = p_img[ix-1 + kk] - bgval;
1428  while ( source > clip && ix <= iend )
1429  {
1430  ifram[1] = ix;
1431  source = p_img[ix-1 + kk] -bgval;
1432  ix ++;
1433  }
1434 
1435  ix = CGN_NINT( xypos[0] );
1436  istr = ifram[2];
1437  source = p_img[ix-1 + kk] - bgval;
1438  while ( source > clip && iy >= istr )
1439  {
1440  ifram[2] = iy;
1441  source = p_img[ix-1 + (*npix *(iy-1))] -bgval;
1442  iy --;
1443  }
1444 
1445  iy = CGN_NINT( xypos[1] );
1446  iend = ifram[3];
1447  source = p_img[ix-1 + kk] - bgval;
1448  while ( source > clip && iy <= iend )
1449  {
1450  ifram[3] = iy;
1451  source = p_img[ix-1 + (*npix *(iy-1))] -bgval;
1452  iy++;
1453  }
1454  nrx = ifram[1] - ifram[0] + 1;
1455  nry = ifram[3] - ifram[2] + 1;
1456  }
1457 
1458  istat = 3; /* iteration failed */
1459 
1460 end_of_iter:
1461  xypos[0] --;
1462  xypos[1] --;
1463  }
1464 
1465 
1466 /* GAUSSIAN centering */
1467 
1468 else
1469  {
1470  register int ii;
1471  int found, ierr, xlim[2], ylim[2], lnew[2];
1472  float lamda, xcent, ycent, xwidth, ywidth;
1473  double chisqr, oldchi, sigma, *krx, *kry, *gfit, *xpos, *yfit,
1474  gpar[MAXPAR];
1475 
1476 
1477  /* construct two marginal distibutions (in pixel coordinates!) */
1478 
1479  lnew[0] = (nry / 4); /* in C notation, from 0 ... */
1480  lnew[1] = nry - (nry / 4) - 1;
1481 
1482 
1483 /* Take care of 1-dim case */
1484 
1485  if (nry == 1)
1486  {
1487  krx = (double *) cpl_calloc(nrx , sizeof(double));
1488  Crhox(p_img,npix,image,lnew,krx);
1489  ierr = Cserch(krx,nrx,IGNORE,xlim,xlim+1,&xcent,&xwidth);
1490 
1491  /* store the data of the fit */
1492 
1493  nval = xlim[1] - xlim[0] + 1;
1494  xpos = (double *) cpl_calloc( nval , sizeof( double ));
1495  yfit = (double *) cpl_calloc( nval , sizeof( double ));
1496  gfit = (double *) cpl_calloc( nval , sizeof( double ));
1497  for (ii=0; ii<nval; ii++)
1498  {
1499  xpos[ii] = xlim[0] + ii;
1500  yfit[ii] = krx[xlim[0] + ii];
1501  }
1502 
1503  /* set parameters for LSQFIT (old FITINTE) */
1504 
1505  lamda = 0.001;
1506  chisqr = GCHIMAX;
1507  gpar[0] = krx[CGN_NINT(xcent)];
1508  gpar[1] = xcent;
1509  gpar[2] = xwidth;
1510  gpar[3] = (krx[xlim[0]] + krx[xlim[1]]) / 2;
1511  (void) cpl_free( (char *) krx );
1512 
1513  it = 0;
1514  found = FALSE;
1515  while ( ! found && it++ < GMXITER )
1516  {
1517  oldchi = chisqr;
1518  ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma);
1519  if ( ierr != 0 )
1520  {
1521  found = NOCONV;
1522  istat = 3;
1523  }
1524  else if ( (oldchi - chisqr)/ chisqr < GCHIFND )
1525  found = TRUE;
1526  }
1527 
1528  /* Is the source still in the image and has it a resonable shape? */
1529 
1530  lamda = 0.0;
1531  ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma);
1532  if ( ierr != 0 )
1533  {
1534  found = NOCONV;
1535  istat = 3;
1536  }
1537  else
1538  {
1539  sumi = (float)(gpar[1] + image[0]);
1540  indx = CGN_NINT(sumi);
1541  if ( indx < 0 || indx >= *npix )
1542  {
1543  found = OUTSIDE;
1544  istat = 2;
1545  }
1546  }
1547  (void) cpl_free( (char *) xpos );
1548  (void) cpl_free( (char *) yfit );
1549  (void) cpl_free( (char *) gfit );
1550 
1551  if ( found == TRUE )
1552  {
1553  xypos[0] = sumi;
1554  xypos[1] = 0;
1555  xysig[0] = (float) gpar[2];
1556  xyerr[0] = (float) sqrt( sigma * chisqr );
1557  indx = CGN_NINT( xypos[0]);
1558  *xyval = p_img[indx];
1559  }
1560  }
1561 
1562 
1563 /* Take care of 2-dim case */
1564 
1565  else
1566  {
1567  krx = (double *) cpl_calloc( nrx , sizeof( double ));
1568  kry = (double *) cpl_calloc( nry , sizeof( double ));
1569 
1570  /* Compute and search X-marginal & Y-marginal */
1571 
1572  Crhox( p_img, npix, image, lnew, krx );
1573  ierr = Cserch( krx, nrx, IGNORE, xlim, xlim+1, &xcent, &xwidth );
1574  lnew[0] = MYMAX( xlim[0], CGN_NINT(xcent - (2 * xwidth)));
1575  lnew[1] = MYMIN( xlim[1], CGN_NINT(xcent + (2 * xwidth)));
1576 
1577  Crhoy( p_img, npix, image, lnew, kry );
1578  ierr = Cserch( kry, nry, IGNORE, ylim, ylim+1, &ycent, &ywidth );
1579  lnew[0] = MYMAX( ylim[0], CGN_NINT(ycent - (2 * ywidth)));
1580  lnew[1] = MYMIN( ylim[1], CGN_NINT(ycent + (2 * ywidth)));
1581 
1582  Crhox( p_img, npix, image, lnew, krx );
1583  ierr = Cserch( krx, nrx, IGNORE, xlim, xlim+1, &xcent, &xwidth );
1584  lnew[0] = MYMAX( xlim[0], CGN_NINT(xcent - (2 * xwidth)));
1585  lnew[1] = MYMIN( xlim[1], CGN_NINT(xcent + (2 * xwidth)));
1586 
1587  Crhoy( p_img, npix, image, lnew, kry );
1588  ierr = Cserch( kry, nry, IGNORE, ylim, ylim+1, &ycent, &ywidth );
1589 
1590  /* fit a gaussian to the source along the X-axis */
1591 
1592  nval = xlim[1] - xlim[0] + 1;
1593  xpos = (double *) cpl_calloc( nval , sizeof( double ));
1594  yfit = (double *) cpl_calloc( nval , sizeof( double ));
1595  gfit = (double *) cpl_calloc( nval , sizeof( double ));
1596  for (ii=0; ii<nval; ii++)
1597  {
1598  xpos[ii] = xlim[0] + ii;
1599  yfit[ii] = krx[xlim[0] + ii];
1600  }
1601 
1602  /* set parameters for LSQFIT */
1603 
1604  lamda = 0.001;
1605  chisqr = GCHIMAX;
1606  gpar[0] = krx[CGN_NINT( xcent )];
1607  gpar[1] = xcent;
1608  gpar[2] = xwidth;
1609  gpar[3] = (krx[xlim[0]] + krx[xlim[1]]) / 2;
1610  (void) cpl_free( ( char *) krx );
1611 
1612  it = 0;
1613  found = FALSE;
1614  while ( ! found && it++ < GMXITER )
1615  {
1616  oldchi = chisqr;
1617  ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma);
1618  if ( ierr != 0 || gpar[2] <= 0.0 )
1619  {
1620  found = NOCONV;
1621  istat = 3;
1622  }
1623  else if ( (oldchi - chisqr)/ chisqr < GCHIFND )
1624  found = TRUE;
1625  }
1626 
1627  /* Is the source still in the image and has it a resonable shape? */
1628 
1629  lamda = 0.0;
1630  ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma);
1631  if ( ierr != 0 )
1632  {
1633  found = NOCONV;
1634  istat = 3;
1635  }
1636  else
1637  {
1638  sumi = (float)(gpar[1] + image[0]);
1639  indx = CGN_NINT(sumi);
1640  if ( indx < 0 || indx >= *npix )
1641  {
1642  found = OUTSIDE;
1643  istat = 2;
1644  }
1645  }
1646 
1647  (void) cpl_free( (char *) xpos );
1648  (void) cpl_free( (char *) yfit );
1649  (void) cpl_free( (char *) gfit );
1650 
1651  if ( found == TRUE )
1652  {
1653  xypos[0] = sumi;
1654  xysig[0] = (float) gpar[2];
1655  xyerr[0] = (float) sqrt( sigma * chisqr );
1656 
1657  /* x-dir o.k. - now fit a gaussian to the source along the Y-axis */
1658 
1659  nval = ylim[1] - ylim[0] + 1;
1660  xpos = (double *) cpl_calloc( nval , sizeof( double ));
1661  yfit = (double *) cpl_calloc( nval , sizeof( double ));
1662  gfit = (double *) cpl_calloc( nval , sizeof( double ));
1663 
1664  for (ii=0; ii<nval; ii++)
1665  {
1666  xpos[ii] = ylim[0] + ii;
1667  yfit[ii] = kry[ylim[0] + ii];
1668  }
1669 
1670  /* set parameters for LSQFIT */
1671 
1672  lamda = 0.001;
1673  chisqr = GCHIMAX;
1674  gpar[0] = kry[CGN_NINT( ycent )];
1675  gpar[1] = ycent;
1676  gpar[2] = ywidth;
1677  gpar[3] = (kry[ylim[0]] + kry[ylim[1]]) / 2;
1678 
1679  it = 0;
1680  found = FALSE;
1681  while ( ! found && it++ < GMXITER )
1682  {
1683  oldchi = chisqr;
1684  ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma );
1685  if ( ierr != 0 || gpar[2] <= 0.0 )
1686  {
1687  found = NOCONV;
1688  istat = 3;
1689  }
1690  else if ( (oldchi - chisqr)/ chisqr < GCHIFND )
1691  found = TRUE;
1692  }
1693 
1694  /* Is the source still in the image and has it a resonable shape? */
1695 
1696  lamda = 0.0;
1697  ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma);
1698  if ( ierr != 0 )
1699  {
1700  found = NOCONV;
1701  istat = 3;
1702  }
1703  else
1704  {
1705  indx = CGN_NINT(xypos[0]);
1706  sumi = (float) (gpar[1] + image[2]);
1707  indy = CGN_NINT(sumi);
1708  if ( indy < 0 || indy >= npix[1] )
1709  {
1710  found = OUTSIDE;
1711  istat = 2;
1712  }
1713  else
1714  indx += (*npix) * indy;
1715  }
1716  (void) cpl_free( (char *) xpos );
1717  (void) cpl_free( (char *) yfit );
1718  (void) cpl_free( (char *) gfit );
1719 
1720  if ( found == TRUE )
1721  {
1722  xypos[1] = sumi;
1723  xysig[1] = (float) gpar[2];
1724  xyerr[1] = (float) sqrt( sigma * chisqr );
1725  *xyval = p_img[indx];
1726  }
1727  }
1728  (void) cpl_free( ( char *) kry );
1729  }
1730  }
1731 
1732 return istat;
1733 }