00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070 #ifdef HAVE_CONFIG_H
00071 # include <config.h>
00072 #endif
00073
00074
00081
00082
00083
00084
00085
00086
00087
00088
00089 #include "uves_physmod_cstacen.h"
00090 #include <cpl.h>
00091 #include <math.h>
00092
00093
00094
00095
00096
00097
00098 #define _POSIX_SOURCE 1
00099
00100
00101
00102 #ifndef PI
00103 #define PI 3.14159265358979325e0
00104 #endif
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114 #define MINVAL 1.0e-37
00115 #define MMXITER 8
00116 #define SMALL 1.0e-20
00117
00118
00119
00120
00121 #define MAXPAR 4
00122 #define IGNORE 2
00123 #define NOCONV -1
00124 #define OUTSIDE -2
00125 #define GMXITER 50
00126 #define GCHIMAX 5.0e+16
00127 #define GCHIFND 0.005
00128
00129
00130 #define MYMIN(a,b) ((a) > (b) ? (b) : (a))
00131 #define MYMAX(a,b) ((b) > (a) ? (b) : (a))
00132
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00185
00186
00187 int
00188 uves_physmod_stacen(float* p_img, int dimx, int dimy, char meth, int* image,
00189 float* xout, float* yout, float* xerr, float* yerr,
00190 float* xsig, float* ysig, float* xyval, int* stat )
00191
00192 {
00193 int npix[2], imap[4];
00194 float xypos[2], xyerr[2], xysig[2];
00195
00196 npix[0] = dimx;
00197 npix[1] = dimy;
00198
00199 imap[0] = image[0] - 1;
00200 imap[1] = image[1] - 1;
00201 imap[2] = image[2] - 1;
00202 imap[3] = image[3] - 1;
00203
00204
00205
00206
00207
00208
00209 *stat = uves_physmod_cstacen(meth, p_img, npix, imap, xypos, xyerr, xysig, xyval );
00210
00211 *xout = xypos[0] + 1;
00212 *yout = xypos[1] + 1;
00213
00214
00215 *xerr = xyerr[0];
00216 *yerr = xyerr[1];
00217 *xsig = xysig[0];
00218 *ysig = xysig[1];
00219
00220
00221
00222
00223
00224 return 0;
00225 }
00226
00232 static int CGN_NINT(float a){
00233 int res=0;
00234
00235 res = (a) > 0 ? floor(a+0.5) : -floor(a-0.5);
00236 return res;
00237 }
00238
00239
00266
00267
00268
00269 static int Ckapsig( float *val, int nval, int iter, float akap,
00270 float *cons, float *rms, int *npts )
00271
00272 {
00273
00274 register int ii=0;
00275 register int it=0;
00276 register int nr=0;
00277
00278 int nr_old=0;
00279
00280 float clip=0;
00281 float dels=0;
00282 float delv=0;
00283 float mean=0;
00284 float msq=0;
00285 float sum=0;
00286 float *vsq=NULL;
00287
00288 if ( nval < 2 ) return (-1);
00289
00290
00291
00292
00293 mean = 0.0;
00294 for (ii=0; ii<nval; ii++) mean += val[ii];
00295 mean /= (float) (nval);
00296 msq = mean * mean;
00297
00298
00299
00300
00301 vsq = (float *) cpl_calloc( nval, sizeof( float ));
00302
00303 dels = 0.0;
00304 for (ii=0; ii<nval; ii++)
00305 {
00306 vsq[ii] = val[ii] * val[ii];
00307 delv = MYMAX( 0.0, vsq[ii] + msq - (2.0 * mean * val[ii]));
00308 dels += delv;
00309 }
00310
00311 *rms = (float) sqrt( MYMAX( MINVAL, dels / (nval-1)));
00312 clip = akap * (*rms);
00313
00314
00315
00316
00317 nr_old = 0;
00318 for (it=0; it<iter; it++)
00319 {
00320 nr = 0;
00321 sum = 0.0;
00322 dels = 0.0;
00323
00324 for ( ii = 0; ii < nval; ii++ )
00325 {
00326 if ( fabs( val[ii] - mean ) < (double) clip )
00327 {
00328 nr++;
00329 delv = MYMAX( 0.0, vsq[ii] + msq - 2.0 * mean * val[ii]);
00330 dels += delv;
00331 sum += val[ii];
00332 }
00333 }
00334
00335 if ( nr <= 2 || nr == nr_old ) goto end_of_it;
00336
00337
00338
00339
00340 nr_old = nr;
00341 *rms = (float) sqrt( MYMAX( MINVAL, dels / (nr-1)));
00342 clip = akap * *rms;
00343 mean = sum / nr_old;
00344 msq = mean * mean;
00345 }
00346
00347 end_of_it:
00348 *cons = mean;
00349 *npts = nr;
00350
00351 cpl_free(vsq);
00352 return (0);
00353 }
00354
00355
00356
00371
00372
00373 static int MATINV( double (*matrix)[MAXPAR], int nfree )
00374 {
00375
00376 register int ii=0;
00377 register int jj=0;
00378 register int kk=0;
00379
00380 int evin=0;
00381 int row=0;
00382 int per[MAXPAR];
00383 double even=0.;
00384 double mjk=0.;
00385 double rowmax=0.;
00386 double hv[MAXPAR];
00387
00388
00389 for ( ii = 0; ii < nfree; ii++ ) per[ii] = ii;
00390
00391 for ( jj = 0; jj < nfree; jj++ )
00392 {
00393 rowmax = fabs( matrix[jj][jj] );
00394 row = jj;
00395 for ( ii = jj + 1; ii < nfree; ii++ )
00396 {
00397 if ( fabs( matrix[ii][jj] ) > rowmax )
00398 {
00399 rowmax = fabs( matrix[ii][jj] );
00400 row = ii;
00401 }
00402 }
00403
00404 if (fabs(matrix[row][jj]) < SMALL)
00405 return (1);
00406
00407 if ( row > jj )
00408 {
00409 for ( kk = 0; kk < nfree; kk++ )
00410 {
00411 even = matrix[jj][kk];
00412 matrix[jj][kk] = matrix[row][kk];
00413 matrix[row][kk] = even;
00414 }
00415 evin = per[jj];
00416 per[jj] = per[row];
00417 per[row] = evin;
00418 }
00419
00420 even = 1.0 / matrix[jj][jj];
00421 for (ii=0; ii<nfree; ii++)
00422 matrix[ii][jj] *= even;
00423 matrix[jj][jj] = even;
00424 for (kk=0; kk<jj; kk++)
00425 {
00426 mjk = matrix[jj][kk];
00427 for ( ii = 0; ii < jj; ii++ )
00428 matrix[ii][kk] -= matrix[ii][jj] * mjk;
00429 for ( ii = jj + 1; ii < nfree; ii++ )
00430 matrix[ii][kk] -= matrix[ii][jj] * mjk;
00431 matrix[jj][kk] = -even * mjk;
00432 }
00433
00434 for ( kk = jj + 1; kk < nfree; kk++ )
00435 {
00436 mjk = matrix[jj][kk];
00437 for ( ii = 0; ii < jj; ii++ )
00438 matrix[ii][kk] -= matrix[ii][jj] * mjk;
00439 for ( ii = jj + 1; ii < nfree; ii++ )
00440 matrix[ii][kk] -= matrix[ii][jj] * mjk;
00441 matrix[jj][kk] = -even * mjk;
00442 }
00443 }
00444
00445 for ( ii = 0; ii < nfree; ii++ )
00446 {
00447 for ( kk = 0; kk < nfree; kk++ )
00448 {
00449 hv[per[kk]] = matrix[ii][kk];
00450 }
00451 for ( kk = 0; kk < nfree; kk++)
00452 {
00453 matrix[ii][kk] = hv[kk];
00454 }
00455 }
00456
00457 return 0;
00458 }
00459
00460
00471
00472
00473 static double ERFCC( double xx )
00474
00475 {
00476
00477 double t=0.;
00478 double z=0.;
00479 double ans=0.;
00480 double zz=0.;
00481 double zozo=0.;
00482
00483
00484 z = fabs( xx );
00485 t = 1.0 / (1.0 + (0.5 * z));
00486
00487
00488
00489
00490
00491
00492
00493
00494
00495
00496
00497
00498
00499
00500
00501
00502
00503 zz = -z * z - 1.26551223 + t * ( 1.00002368 + t *
00504 ( 0.37409196 + t * ( 0.09678418 + t *
00505 (-0.18628806 + t * ( 0.27886807 + t *
00506 (-1.13520398 + t * ( 1.48851587 + t *
00507 (-0.82215223 + t * 0.17087277 ))))))));
00508
00509
00510 if (zz < -500.0)
00511 zozo = 0.0;
00512 else
00513 zozo = exp(zz);
00514
00515 ans = t * zozo;
00516
00517
00518 return (xx >= 0.0 ? ans : 2.0 - ans);
00519 }
00520
00521
00531
00532
00533 static double GAUSFU( double xx, double *gpar )
00534 {
00535
00536 double rc=0.;
00537 double dd=0.;
00538
00539 static int init = TRUE;
00540 static double sqrt_2;
00541 static double rc1;
00542
00543
00544 if ( init )
00545 {
00546 sqrt_2 = sqrt( 2.0 );
00547 rc1 = sqrt(PI)/sqrt_2;
00548 init = FALSE;
00549 }
00550
00551 rc = 1.0 / (sqrt_2 * gpar[2]);
00552 dd = xx - gpar[1];
00553 dd = ERFCC(rc * (dd - 0.5)) - ERFCC(rc * (dd + 0.5));
00554 return ( gpar[3] + rc1 * gpar[0] * gpar[2] * dd );
00555 }
00556
00557
00569
00570
00571 static void GAUSDE( double xdat, double *gpar, double *deriv )
00572 {
00573 double temp=0.;
00574 double tempp=0.;
00575 double x1=0.;
00576 double x2=0.;
00577 double zz=0.;
00578 double zx=0.;
00579 double dv1=0.;
00580 double dv2=0.;
00581 static double sqrt_2;
00582
00583 register int jj=0;
00584 static int init = TRUE;
00585
00586
00587 if ( init )
00588 {
00589 sqrt_2 = sqrt( 2.0 );
00590 init = FALSE;
00591 }
00592
00593 temp = sqrt_2 * gpar[2];
00594 tempp = xdat - gpar[1];
00595 x1 = (tempp - 0.5) / temp;
00596 x2 = (tempp + 0.5) / temp;
00597 zz = tempp / gpar[2] ;
00598
00599 if ( ((zz * zz) - 50.0) < 0.0 )
00600 {
00601 deriv[0] = (GAUSFU( xdat, gpar ) - gpar[3]) / gpar[0];
00602
00603 zx = (-x1) * x1;
00604 if ( zx < -200.0)
00605 dv1 = 0.0;
00606 else
00607 dv1 = exp(zx);
00608 zx = (-x2) * x2;
00609 if ( zx < -200.0)
00610 dv2 = dv1;
00611 else
00612 dv2 = dv1 - exp(zx);
00613 deriv[1] = gpar[0] * dv2;
00614
00615
00616
00617
00618
00619 deriv[2] = deriv[1] * zz;
00620 }
00621 else
00622 for (jj=0; jj<3; jj++) deriv[jj] = 0.0;
00623
00624 deriv[3] = 1.0;
00625 }
00626
00627
00643
00644
00645 static float FCHIS(double *data,int ndim,int nfree,int mode,double *dfit) {
00646
00647 register int ii=0;
00648
00649 double diff=0.;
00650 double weight=0.;
00651 double chisq=0.;
00652
00653
00654 if ( nfree > 0 )
00655 {
00656 chisq = 0.0;
00657
00658 for (ii=0; ii<ndim; ii++)
00659 {
00660 if ( mode < 0 )
00661 {
00662 if ( *data < 0 )
00663 weight = -1. / *data;
00664 else if ( *data == 0 )
00665 weight = 1.0;
00666 else
00667 weight = 1. / *data;
00668 }
00669 else
00670 weight = 1.0;
00671
00672 diff = (*data) - (*dfit);
00673 data++; dfit++;
00674 chisq += weight * diff * diff;
00675 }
00676 return (chisq / nfree);
00677 }
00678
00679 else
00680 return 0.0;
00681 }
00682
00683
00711 static int LSQFIT( double *xdat, double *data, int ndim,
00712 double *gpar, float *lamda, double *dfit,
00713 double *chisqr, double *sigma )
00714
00715 {
00716
00717 register int icnt=0;
00718 register int ii=0;
00719 register int jj=0;
00720 register int kk=0;
00721
00722 int nfree=0;
00723
00724 double chisq1=0;
00725 double b[MAXPAR], beta[MAXPAR], deriv[MAXPAR],
00726 array[MAXPAR][MAXPAR], alpha[MAXPAR][MAXPAR];
00727
00728
00729 nfree = ndim - MAXPAR;
00730 *sigma = 0.0;
00731 if ( nfree < 1 || fabs( (double) *gpar ) < SMALL ) return (1);
00732
00733
00734
00735
00736 for (ii=0; ii<MAXPAR; ii++)
00737 {
00738 beta[ii] = 0.0;
00739 for (jj=0; jj<=ii; jj++) alpha[ii][jj] = 0.0;
00740 }
00741
00742 for (ii=0; ii<ndim; ii++)
00743 {
00744 GAUSDE( xdat[ii], gpar, deriv );
00745
00746 for (jj=0; jj<MAXPAR; jj++)
00747 {
00748 beta[jj] += (data[ii] - GAUSFU( xdat[ii], gpar )) * deriv[jj];
00749 for (kk=0; kk<=jj; kk++)
00750 alpha[jj][kk] += deriv[jj] * deriv[kk];
00751 }
00752 }
00753
00754 for (ii=0; ii<MAXPAR; ii++)
00755 {
00756 for (jj=0; jj<=ii; jj++)
00757 alpha[jj][ii] = alpha[ii][jj];
00758 }
00759
00760
00761
00762
00763 if ( *lamda < SMALL)
00764 {
00765 if (MATINV(alpha,MAXPAR) == 1) return (2);
00766
00767 *sigma = MYMAX( 0.0, alpha[1][1] );
00768 }
00769
00770 else
00771 {
00772 for (ii=0; ii<ndim; ii++)
00773 dfit[ii] = GAUSFU( xdat[ii], gpar );
00774
00775 chisq1 = FCHIS( data, ndim, nfree, 0, dfit );
00776
00777 icnt = 0;
00778 loop:
00779 for ( jj = 0; jj < MAXPAR; jj++ )
00780 {
00781 for ( kk = 0; kk < MAXPAR; kk++ )
00782 {
00783 if (fabs( alpha[jj][jj] ) < 1.e-15 || fabs( alpha[kk][kk] ) < 1.e-15)
00784 return 2;
00785 array[jj][kk] = alpha[jj][kk] /
00786 sqrt( alpha[jj][jj] * alpha[kk][kk] ) ;
00787 }
00788 array[jj][jj] = 1.0 + *lamda;
00789 }
00790
00791 (void) MATINV( array, MAXPAR );
00792
00793 for ( jj = 0; jj < MAXPAR; jj++ )
00794 {
00795 b[jj] = gpar[jj] ;
00796 for ( kk = 0; kk < MAXPAR ; kk++ )
00797 {
00798 b[jj] += beta[kk] * array[jj][kk] /
00799 sqrt( alpha[jj][jj] * alpha[kk][kk] );
00800 }
00801 }
00802
00803
00804
00805 for (ii=0; ii<ndim; ii++)
00806 dfit[ii] = GAUSFU( xdat[ii], b );
00807
00808 *chisqr = FCHIS( data, ndim, nfree, 0, dfit );
00809
00810 if ( chisq1 - *chisqr < 0.0 )
00811 {
00812 if (++icnt < 60)
00813 {
00814 *lamda *= 10;
00815 goto loop;
00816 }
00817 else
00818 return (2);
00819 }
00820
00821 for (jj=0; jj<MAXPAR; jj++) gpar[jj] = b[jj];
00822 *lamda /= 10.0;
00823 }
00824
00825 return 0;
00826 }
00827
00828
00829
00849 static void Crhox( float *p_img, int *npix, int *image,
00850 int *lnew, double *krx )
00851
00852 {
00853 register int nxdim=0;
00854 register int ix=0;
00855 register int iy=0;
00856
00857 int nrx=0;
00858 int nry=0;
00859
00860 double sum=0;
00861
00862
00863
00864
00865
00866
00867
00868
00869
00870
00871
00872
00873
00874
00875
00876
00877
00878
00879
00880
00881
00882
00883
00884
00885 nrx = image[1] - image[0] + 1;
00886 nry = lnew[1] - lnew[0] + 1;
00887 nxdim = *npix;
00888 p_img += nxdim * (image[2] + lnew[0]);
00889
00890 for (ix=0; ix<nrx; ix++)
00891 {
00892 sum = 0.0;
00893 for (iy=0; iy<nry*nxdim; iy+=nxdim) sum += p_img[iy];
00894 p_img++;
00895 *krx++ = sum;
00896 }
00897 }
00898
00899
00900
00901
00921 static void Crhoy( float *p_img, int *npix, int *image,
00922 int *lnew, double *kry )
00923
00924 {
00925 register int nxdim=0;
00926 register int ix=0;
00927 register int iy=0;
00928
00929 int nrx=0;
00930 int nry=0;
00931
00932 double sum=0;
00933
00934
00935
00936
00937
00938
00939
00940
00941
00942
00943
00944
00945
00946
00947
00948
00949
00950
00951
00952
00953
00954
00955
00956 nrx = lnew[1] - lnew[0] + 1;
00957 nry = image[3] - image[2] + 1;
00958 nxdim = *npix;
00959 p_img += (nxdim * image[2]) + (image[0] + lnew[0]);
00960
00961 for (iy=0; iy<nry; iy++)
00962 {
00963 sum = 0.0;
00964 for (ix=0; ix<nrx; ix++) sum += p_img[ix];
00965 p_img += nxdim;
00966 *kry++ = sum;
00967 }
00968 }
00969
00970
00994 static int Cserch( double *marg, int ndim, int ign,
00995 int *lmin, int *lmax, float *s_cent, float *s_width )
00996
00997 {
00998
00999 register int ii=0;
01000
01001 int ql=0;
01002 int ibgn=0;
01003 int icrowd=0;
01004 int iend=0;
01005 int imax=0;
01006 int imin=0;
01007 int indx=0;
01008
01009 double dxk=0.;
01010 double diff=0.;
01011 double drmn=0.;
01012 double drmx=0.;
01013 double sum=0.;
01014 double *work=NULL;
01015
01016
01017 ibgn = ign;
01018 iend = ndim - ign -1;
01019
01020
01021
01022 work = (double *) cpl_calloc( ndim , sizeof( double ));
01023
01024
01025
01026
01027 imin = imax = 0;
01028 drmn = drmx = 0.0;
01029 for (ii = ibgn; ii < iend; ii++ )
01030 {
01031 diff = marg[ii+1] - marg[ii-1];
01032 work[ii] = marg[ii+2] - marg[ii-2] + (2 * diff);
01033 if ( work[ii] >= drmx )
01034 {
01035 drmx = work[ii];
01036 imax = ii;
01037 }
01038 if (
01039 work[ii] <= drmn )
01040 { drmn = work[ii];
01041 imin = ii;
01042 }
01043 }
01044
01045
01046
01047
01048 icrowd = 0;
01049 if (imin <= imax )
01050 {
01051 if ( ndim - imax > imin )
01052 {
01053 icrowd = -1;
01054 drmn = drmx;
01055 for ( ii = imax+1; ii < iend; ii++ )
01056 {
01057 if ( work[ii] < drmn )
01058 {
01059 drmn = work[ii];
01060 imin = ii;
01061 }
01062 }
01063 }
01064 else
01065 {
01066 icrowd = 1;
01067 drmx = drmn;
01068 for ( ii = ibgn; ii < imin; ii++ )
01069 {
01070 if ( work[ii] >= drmx )
01071 {
01072 drmx = work[ii];
01073 imax = ii;
01074 }
01075 }
01076 }
01077 }
01078
01079
01080
01081
01082 *s_cent = ((float)(imax + imin)) * 0.5;
01083 *s_width = imin - imax;
01084
01085 sum = 0.0;
01086 for ( ii = imax; ii <= imin; ii++ ) sum += work[ii];
01087
01088 diff = drmx - drmn;
01089 if ( fabs(diff) > SMALL)
01090 {
01091 dxk = sum * *s_width / ( (*s_width+1.0)*diff );
01092 *s_cent += dxk;
01093 }
01094 *s_width /= 2;
01095 indx = CGN_NINT(*s_cent);
01096 if (indx < 0)
01097 {
01098 *s_cent = 0.0;
01099 indx = 0;
01100 }
01101 else if (indx >= ndim)
01102 {
01103 *s_cent = (float)(ndim-1);
01104 indx = ndim-1;
01105 }
01106
01107
01108
01109
01110 ql = indx - 2;
01111 *lmin = 0;
01112 if (ql < 2) goto next_step;
01113
01114 low_loop:
01115 ql --;
01116 if (ql <= 0 ) goto next_step;
01117 if (ql == 1) goto lo5;
01118 if (ql == 2) goto lo4;
01119 if (ql == 3) goto lo3;
01120
01121 if (marg[ql] > marg[ql-4]) goto low_loop;
01122 lo3:
01123 if (marg[ql] > marg[ql-3]) goto low_loop;
01124 lo4:
01125 if (marg[ql] > marg[ql-2]) goto low_loop;
01126 lo5:
01127 if (marg[ql] > marg[ql-1]) goto low_loop;
01128
01129 *lmin = ql + 1;
01130
01131
01132
01133
01134 next_step:
01135 ql = indx + 2;
01136 *lmax = ndim - 1;
01137 ii = ndim - ql;
01138 if (ii < 3) goto end_of_it;
01139
01140 hi_loop:
01141 ql ++ ;
01142 ii = ndim - ql;
01143 if (ii == 1 ) goto end_of_it;
01144 if (ii == 2) goto hi5;
01145 if (ii == 3) goto hi4;
01146 if (ii == 4) goto hi3;
01147
01148 if (marg[ql] > marg[ql+4]) goto hi_loop;
01149 hi3:
01150 if (marg[ql] > marg[ql+3]) goto hi_loop;
01151 hi4:
01152 if (marg[ql] > marg[ql+2]) goto hi_loop;
01153 hi5:
01154 if (marg[ql] > marg[ql+1]) goto hi_loop;
01155 *lmax = ql - 1;
01156
01157 end_of_it:
01158 (void) cpl_free( (char *) work );
01159 return (icrowd);
01160
01161 }
01162
01163
01201
01202
01203
01204
01205
01206
01207
01208 int uves_physmod_cstacen(char meth, float* p_img, int* npix, int* image,
01209 float* xypos, float* xyerr, float* xysig, float* xyval )
01210
01211 {
01212
01213 register int it=0;
01214 register int ix=0;
01215 register int iy=0;
01216
01217 int bgnr=0;
01218 int indx=0;
01219 int indy=0;
01220 int istat=0;
01221 int nrx=0;
01222 int nry=0;
01223 int nval=0;
01224 int ifram[4];
01225
01226 float bgval=0.;
01227 float clip=0.;
01228 float rms=0.;
01229 float xmom=0.;
01230 float ymom=0.;
01231 float source=0.;
01232 float sumi=0.;
01233 float xold=0.;
01234 float yold=0.;
01235
01236 float *p_bgn=NULL;
01237 float *p_edge=NULL;
01238
01239
01240 istat = 0;
01241 p_bgn = p_img;
01242 for (ix=0; ix<4; ix++)
01243 ifram[ix] = image[ix] + 1;
01244 nrx = ifram[1] - ifram[0] + 1;
01245 nry = ifram[3] - ifram[2] + 1;
01246
01247 xypos[0] = (ifram[0] + ifram[1]) * 0.5;
01248 xypos[1] = (ifram[2] + ifram[3]) * 0.5;
01249 xyerr[0] = xyerr[1] = 0.0;
01250 xysig[0] = xysig[1] = 0.0;
01251 *xyval = 0.0;
01252
01253
01254
01255 if ( meth != 'G' && meth != 'g' )
01256 {
01257 int kk, istr, iend;
01258
01259 xold = yold = -1.0;
01260
01261 p_img += (ifram[0] - 1) + (npix[0] * (ifram[2] - 1));
01262
01263
01264
01265 if (nry > 1)
01266 {
01267 nval = (2 * nrx) + (2 * (nry-2));
01268 p_edge = (float *) cpl_calloc( nval , sizeof( float ));
01269
01270 for (ix=0; ix<nrx;ix++)
01271 *p_edge++ = p_img[ix];
01272
01273 p_img += *npix;
01274 for (iy=0; iy<(nry-2); iy++)
01275 {
01276 *p_edge++ = p_img[0];
01277 *p_edge++ = p_img[nrx - 1];
01278 p_img += npix[0];
01279 }
01280 for (ix=0; ix<nrx; ix++) *p_edge++ = p_img[ix];
01281 }
01282 else
01283 {
01284 nval = nrx;
01285 p_edge = (float *) cpl_calloc( nval , sizeof( float ));
01286
01287 for (ix=0; ix<nrx;ix++)
01288 *p_edge++ = p_img[ix];
01289 }
01290
01291 p_img = p_bgn;
01292
01293 p_edge -= nval;
01294 (void) Ckapsig( p_edge, nval, 5, 2.0, &bgval, &rms, &bgnr );
01295 (void) cpl_free( (char *) p_edge );
01296
01297
01298
01299
01300 clip = 3.0 * rms;
01301
01302 for (it=0; it<MMXITER; it++)
01303 {
01304 sumi = xmom = ymom = 0.0;
01305 p_img += ifram[0] - 1 + (npix[0] * (ifram[2] - 1));
01306 for (iy=0; iy<nry; iy++)
01307 {
01308 for (ix=0; ix<nrx; ix++)
01309 {
01310 if ( (source = p_img[ix] - bgval) > clip )
01311 {
01312 sumi += source;
01313 xmom += source * (ifram[0] + ix);
01314 ymom += source * (ifram[2] + iy);
01315 }
01316 }
01317 p_img += npix[0];
01318 }
01319 p_img = p_bgn;
01320
01321 if ((nrx < 3) || (nry < 3))
01322 {
01323 xysig[0] = nrx;
01324 xysig[1] = nry;
01325 istat = 1;
01326 if ( sumi > 0.0 )
01327 {
01328 xypos[0] = xmom / sumi;
01329 xypos[1] = ymom / sumi;
01330 }
01331 else
01332 {
01333 istat = 2;
01334 xypos[0] = (ifram[0] + ifram[1]) * 0.5;
01335 xypos[1] = (ifram[2] + ifram[3]) * 0.5;
01336 }
01337 indx = CGN_NINT(xypos[0]-1);
01338 indy = CGN_NINT(xypos[1]-1);
01339 *xyval = p_img[indx + ((*npix) * indy)];
01340 goto end_of_iter;
01341 }
01342
01343 if (sumi > 0.0)
01344 {
01345 xypos[0] = xmom / sumi;
01346 xypos[1] = ymom / sumi;
01347 xysig[0] = nrx;
01348 xysig[1] = nry;
01349
01350 if ( xold == xypos[0] && yold == xypos[1] )
01351 {
01352 int nr = 0;
01353 double xdif, ydif, xrms, yrms;
01354
01355 xrms = yrms = sumi = 0.0;
01356 p_img += ifram[0] - 1 + (npix[0] * (ifram[2] - 1));
01357 for (iy=0; iy<nry; iy++ )
01358 {
01359 for (ix=0; ix<nrx; ix++)
01360 {
01361 if ( (source = p_img[ix] - bgval) > clip )
01362 {
01363 xdif = (ifram[0] + ix) - xypos[0];
01364 ydif = (ifram[2] + iy) - xypos[1];
01365 xrms += fabs( source * xdif *xdif );
01366 yrms += fabs( source * ydif *ydif );
01367 sumi += source;
01368 nr++;
01369 }
01370 }
01371 p_img += npix[0];
01372 }
01373 p_img = p_bgn;
01374
01375 indx = CGN_NINT(xypos[0]-1) + (npix[0] * CGN_NINT(xypos[1]-1));
01376 *xyval = p_img[indx];
01377 xysig[0] = (float) sqrt(xrms /(sumi+ *xyval - bgval));
01378 xysig[1] = (float) sqrt(yrms /(sumi+ *xyval - bgval));
01379 xyerr[0] = (float) (xysig[0] / sqrt( (double) (nr - 1)));
01380 xyerr[1] = (float) (xysig[1] / sqrt( (double) (nr - 1)));
01381 goto end_of_iter;
01382 }
01383
01384
01385 xold = xypos[0];
01386 yold = xypos[1];
01387 }
01388 else
01389 {
01390 istat = 2;
01391 xypos[0] = (ifram[0] + ifram[1]) * 0.5;
01392 xypos[1] = (ifram[2] + ifram[3]) * 0.5;
01393 indx = CGN_NINT(xypos[0]-1);
01394 indy = CGN_NINT(xypos[1]-1);
01395 *xyval = p_img[indx + ((*npix) * indy)];
01396 goto end_of_iter;
01397 }
01398
01399
01400
01401
01402 indx = CGN_NINT(xypos[0]-1) + (npix[0] * CGN_NINT(xypos[1]-1));
01403 if ( (*xyval = p_img[indx] - bgval) <= clip )
01404 {
01405 xysig[0] = xysig[1] = 0.0;
01406 istat = 1;
01407 goto end_of_iter;
01408 }
01409
01410
01411
01412
01413 ix = CGN_NINT( xypos[0] );
01414 iy = CGN_NINT( xypos[1] );
01415 kk = npix[0] * (iy - 1);
01416 istr = ifram[0];
01417 source = p_img[ix-1 + kk] - bgval;
01418 while ( source > clip && ix >= istr )
01419 {
01420 ifram[0] = ix;
01421 source = p_img[ix-1 + kk] - bgval;
01422 ix --;
01423 }
01424
01425 ix = CGN_NINT( xypos[0] );
01426 iend = ifram[1];
01427 source = p_img[ix-1 + kk] - bgval;
01428 while ( source > clip && ix <= iend )
01429 {
01430 ifram[1] = ix;
01431 source = p_img[ix-1 + kk] -bgval;
01432 ix ++;
01433 }
01434
01435 ix = CGN_NINT( xypos[0] );
01436 istr = ifram[2];
01437 source = p_img[ix-1 + kk] - bgval;
01438 while ( source > clip && iy >= istr )
01439 {
01440 ifram[2] = iy;
01441 source = p_img[ix-1 + (*npix *(iy-1))] -bgval;
01442 iy --;
01443 }
01444
01445 iy = CGN_NINT( xypos[1] );
01446 iend = ifram[3];
01447 source = p_img[ix-1 + kk] - bgval;
01448 while ( source > clip && iy <= iend )
01449 {
01450 ifram[3] = iy;
01451 source = p_img[ix-1 + (*npix *(iy-1))] -bgval;
01452 iy++;
01453 }
01454 nrx = ifram[1] - ifram[0] + 1;
01455 nry = ifram[3] - ifram[2] + 1;
01456 }
01457
01458 istat = 3;
01459
01460 end_of_iter:
01461 xypos[0] --;
01462 xypos[1] --;
01463 }
01464
01465
01466
01467
01468 else
01469 {
01470 register int ii;
01471 int found, ierr, xlim[2], ylim[2], lnew[2];
01472 float lamda, xcent, ycent, xwidth, ywidth;
01473 double chisqr, oldchi, sigma, *krx, *kry, *gfit, *xpos, *yfit,
01474 gpar[MAXPAR];
01475
01476
01477
01478
01479 lnew[0] = (nry / 4);
01480 lnew[1] = nry - (nry / 4) - 1;
01481
01482
01483
01484
01485 if (nry == 1)
01486 {
01487 krx = (double *) cpl_calloc(nrx , sizeof(double));
01488 Crhox(p_img,npix,image,lnew,krx);
01489 ierr = Cserch(krx,nrx,IGNORE,xlim,xlim+1,&xcent,&xwidth);
01490
01491
01492
01493 nval = xlim[1] - xlim[0] + 1;
01494 xpos = (double *) cpl_calloc( nval , sizeof( double ));
01495 yfit = (double *) cpl_calloc( nval , sizeof( double ));
01496 gfit = (double *) cpl_calloc( nval , sizeof( double ));
01497 for (ii=0; ii<nval; ii++)
01498 {
01499 xpos[ii] = xlim[0] + ii;
01500 yfit[ii] = krx[xlim[0] + ii];
01501 }
01502
01503
01504
01505 lamda = 0.001;
01506 chisqr = GCHIMAX;
01507 gpar[0] = krx[CGN_NINT(xcent)];
01508 gpar[1] = xcent;
01509 gpar[2] = xwidth;
01510 gpar[3] = (krx[xlim[0]] + krx[xlim[1]]) / 2;
01511 (void) cpl_free( (char *) krx );
01512
01513 it = 0;
01514 found = FALSE;
01515 while ( ! found && it++ < GMXITER )
01516 {
01517 oldchi = chisqr;
01518 ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma);
01519 if ( ierr != 0 )
01520 {
01521 found = NOCONV;
01522 istat = 3;
01523 }
01524 else if ( (oldchi - chisqr)/ chisqr < GCHIFND )
01525 found = TRUE;
01526 }
01527
01528
01529
01530 lamda = 0.0;
01531 ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma);
01532 if ( ierr != 0 )
01533 {
01534 found = NOCONV;
01535 istat = 3;
01536 }
01537 else
01538 {
01539 sumi = (float)(gpar[1] + image[0]);
01540 indx = CGN_NINT(sumi);
01541 if ( indx < 0 || indx >= *npix )
01542 {
01543 found = OUTSIDE;
01544 istat = 2;
01545 }
01546 }
01547 (void) cpl_free( (char *) xpos );
01548 (void) cpl_free( (char *) yfit );
01549 (void) cpl_free( (char *) gfit );
01550
01551 if ( found == TRUE )
01552 {
01553 xypos[0] = sumi;
01554 xypos[1] = 0;
01555 xysig[0] = (float) gpar[2];
01556 xyerr[0] = (float) sqrt( sigma * chisqr );
01557 indx = CGN_NINT( xypos[0]);
01558 *xyval = p_img[indx];
01559 }
01560 }
01561
01562
01563
01564
01565 else
01566 {
01567 krx = (double *) cpl_calloc( nrx , sizeof( double ));
01568 kry = (double *) cpl_calloc( nry , sizeof( double ));
01569
01570
01571
01572 Crhox( p_img, npix, image, lnew, krx );
01573 ierr = Cserch( krx, nrx, IGNORE, xlim, xlim+1, &xcent, &xwidth );
01574 lnew[0] = MYMAX( xlim[0], CGN_NINT(xcent - (2 * xwidth)));
01575 lnew[1] = MYMIN( xlim[1], CGN_NINT(xcent + (2 * xwidth)));
01576
01577 Crhoy( p_img, npix, image, lnew, kry );
01578 ierr = Cserch( kry, nry, IGNORE, ylim, ylim+1, &ycent, &ywidth );
01579 lnew[0] = MYMAX( ylim[0], CGN_NINT(ycent - (2 * ywidth)));
01580 lnew[1] = MYMIN( ylim[1], CGN_NINT(ycent + (2 * ywidth)));
01581
01582 Crhox( p_img, npix, image, lnew, krx );
01583 ierr = Cserch( krx, nrx, IGNORE, xlim, xlim+1, &xcent, &xwidth );
01584 lnew[0] = MYMAX( xlim[0], CGN_NINT(xcent - (2 * xwidth)));
01585 lnew[1] = MYMIN( xlim[1], CGN_NINT(xcent + (2 * xwidth)));
01586
01587 Crhoy( p_img, npix, image, lnew, kry );
01588 ierr = Cserch( kry, nry, IGNORE, ylim, ylim+1, &ycent, &ywidth );
01589
01590
01591
01592 nval = xlim[1] - xlim[0] + 1;
01593 xpos = (double *) cpl_calloc( nval , sizeof( double ));
01594 yfit = (double *) cpl_calloc( nval , sizeof( double ));
01595 gfit = (double *) cpl_calloc( nval , sizeof( double ));
01596 for (ii=0; ii<nval; ii++)
01597 {
01598 xpos[ii] = xlim[0] + ii;
01599 yfit[ii] = krx[xlim[0] + ii];
01600 }
01601
01602
01603
01604 lamda = 0.001;
01605 chisqr = GCHIMAX;
01606 gpar[0] = krx[CGN_NINT( xcent )];
01607 gpar[1] = xcent;
01608 gpar[2] = xwidth;
01609 gpar[3] = (krx[xlim[0]] + krx[xlim[1]]) / 2;
01610 (void) cpl_free( ( char *) krx );
01611
01612 it = 0;
01613 found = FALSE;
01614 while ( ! found && it++ < GMXITER )
01615 {
01616 oldchi = chisqr;
01617 ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma);
01618 if ( ierr != 0 || gpar[2] <= 0.0 )
01619 {
01620 found = NOCONV;
01621 istat = 3;
01622 }
01623 else if ( (oldchi - chisqr)/ chisqr < GCHIFND )
01624 found = TRUE;
01625 }
01626
01627
01628
01629 lamda = 0.0;
01630 ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma);
01631 if ( ierr != 0 )
01632 {
01633 found = NOCONV;
01634 istat = 3;
01635 }
01636 else
01637 {
01638 sumi = (float)(gpar[1] + image[0]);
01639 indx = CGN_NINT(sumi);
01640 if ( indx < 0 || indx >= *npix )
01641 {
01642 found = OUTSIDE;
01643 istat = 2;
01644 }
01645 }
01646
01647 (void) cpl_free( (char *) xpos );
01648 (void) cpl_free( (char *) yfit );
01649 (void) cpl_free( (char *) gfit );
01650
01651 if ( found == TRUE )
01652 {
01653 xypos[0] = sumi;
01654 xysig[0] = (float) gpar[2];
01655 xyerr[0] = (float) sqrt( sigma * chisqr );
01656
01657
01658
01659 nval = ylim[1] - ylim[0] + 1;
01660 xpos = (double *) cpl_calloc( nval , sizeof( double ));
01661 yfit = (double *) cpl_calloc( nval , sizeof( double ));
01662 gfit = (double *) cpl_calloc( nval , sizeof( double ));
01663
01664 for (ii=0; ii<nval; ii++)
01665 {
01666 xpos[ii] = ylim[0] + ii;
01667 yfit[ii] = kry[ylim[0] + ii];
01668 }
01669
01670
01671
01672 lamda = 0.001;
01673 chisqr = GCHIMAX;
01674 gpar[0] = kry[CGN_NINT( ycent )];
01675 gpar[1] = ycent;
01676 gpar[2] = ywidth;
01677 gpar[3] = (kry[ylim[0]] + kry[ylim[1]]) / 2;
01678
01679 it = 0;
01680 found = FALSE;
01681 while ( ! found && it++ < GMXITER )
01682 {
01683 oldchi = chisqr;
01684 ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma );
01685 if ( ierr != 0 || gpar[2] <= 0.0 )
01686 {
01687 found = NOCONV;
01688 istat = 3;
01689 }
01690 else if ( (oldchi - chisqr)/ chisqr < GCHIFND )
01691 found = TRUE;
01692 }
01693
01694
01695
01696 lamda = 0.0;
01697 ierr = LSQFIT(xpos,yfit,nval,gpar,&lamda,gfit,&chisqr,&sigma);
01698 if ( ierr != 0 )
01699 {
01700 found = NOCONV;
01701 istat = 3;
01702 }
01703 else
01704 {
01705 indx = CGN_NINT(xypos[0]);
01706 sumi = (float) (gpar[1] + image[2]);
01707 indy = CGN_NINT(sumi);
01708 if ( indy < 0 || indy >= npix[1] )
01709 {
01710 found = OUTSIDE;
01711 istat = 2;
01712 }
01713 else
01714 indx += (*npix) * indy;
01715 }
01716 (void) cpl_free( (char *) xpos );
01717 (void) cpl_free( (char *) yfit );
01718 (void) cpl_free( (char *) gfit );
01719
01720 if ( found == TRUE )
01721 {
01722 xypos[1] = sumi;
01723 xysig[1] = (float) gpar[2];
01724 xyerr[1] = (float) sqrt( sigma * chisqr );
01725 *xyval = p_img[indx];
01726 }
01727 }
01728 (void) cpl_free( ( char *) kry );
01729 }
01730 }
01731
01732 return istat;
01733 }