imcore_overlp.c

00001 /* $Id: imcore_overlp.c,v 1.9 2010/10/05 09:11:33 jim Exp $
00002  *
00003  * This file is part of the VIRCAM Pipeline
00004  * Copyright (C) 2005 Cambridge Astronomy Survey Unit
00005  *
00006  * This program is free software; you can redistribute it and/or modify
00007  * it under the terms of the GNU General Public License as published by
00008  * the Free Software Foundation; either version 2 of the License, or
00009  * (at your option) any later version.
00010  *
00011  * This program is distributed in the hope that it will be useful,
00012  * but WITHOUT ANY WARRANTY; without even the implied warranty of
00013  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00014  * GNU General Public License for more details.
00015  *
00016  * You should have received a copy of the GNU General Public License
00017  * along with this program; if not, write to the Free Software
00018  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00019  */
00020 
00021 /*
00022  * $Author: jim $
00023  * $Date: 2010/10/05 09:11:33 $
00024  * $Revision: 1.9 $
00025  * $Name: v1-1-0 $
00026  */
00027 
00028 #include <stdio.h>
00029 #include <stdlib.h>
00030 #include <math.h>
00031 
00032 #include <cpl.h>
00033 
00034 #include "ap.h"
00035 #include "util.h"
00036 #include "imcore.h"
00037 #include "floatmath.h"
00038 
00039 #define NITER 6
00040 
00041 static void moments_thr(ap_t *, float [], int []);
00042 static void sort_on_zsm_rev(int, plstruct *);
00043 static void update_ov(float [], float, float, float, float);
00044 static void check_term(ap_t *, int *, float [IMNUM][NPAR+1], int [IMNUM][2], 
00045                        int *);
00046 
00047 static float oldthr;
00048 static float curthr;
00049 static float nexthr;
00050 static float lasthr;
00051 static float xbar_start;
00052 static float ybar_start;
00053 
00056 /*---------------------------------------------------------------------------*/
00096 /*---------------------------------------------------------------------------*/
00097 
00098 extern void overlp(ap_t *ap, float parm[IMNUM][NPAR], int *nbit, float xbar,
00099                    float ybar, float total, int npix, float tmax) {
00100     plstruct *pl;
00101     int npl,ipix,ipixo2,npl2,nbitprev,nobj,toomany,i,isnew,k,kk,j,ibitx[IMNUM];
00102     int ibity[IMNUM],ibitl[IMNUM],iwas,iupdate[IMNUM],npl3,iter,lastone,ic,jj;
00103     int ii,conv,ipks[IMNUM][2];
00104     float fconst,offset,tmul,smul,xintmn,itmaxlim,algthr,radmax,xb,yb,radius2;
00105     float results[IMNUM][NPAR+1],distmax,dx,dy,parmnew[IMNUM][NPAR],sumint;
00106     float xlevol,radold,slope,xx,xlevel,radius,xdat[NAREAL+1],xcor[NAREAL+1];
00107     float dlbydr,wt,dist,xeff,polycf[3],ttt,radthr,delb,deli,ratio;
00108     float bitx[IMNUM],bitl[IMNUM],sxx,syy;
00109     ap_t ap2;
00110 
00111     /* Initialise a few variables */
00112 
00113     pl = ap->plarray;
00114     npl = ap->npl_pix;
00115     ipix = ap->ipnop;
00116     oldthr = ap->thresh;
00117     fconst = ap->fconst;
00118     offset = ap->areal_offset;
00119     xbar_start = xbar;
00120     ybar_start = ybar;
00121 
00122     /* Initialise some constants that you might need later */
00123 
00124     tmul = 1.2589678;             /* 1/4 mag deblending contour increment */
00125     smul = 2.5;                   /* starting contour increment */
00126     ipixo2 = MAX(2,(ipix + 1)/2); /* ipix is the minimum image pixel size */
00127     xintmn = oldthr*ipixo2;       /* minimum intensity for fragments */
00128     itmaxlim = 0.9*tmax;          /* upper limit deblending 90% of peak */
00129     lasthr = itmaxlim;
00130     algthr = logf(oldthr);        /* Convenient bit of shorthand */
00131     radmax = sqrtf(((float)npix)/CPL_MATH_PI); /* Isophotal size of input data array */
00132 
00133     /* Get a maximum of IDBLIM points brighter than the new detection threshold
00134        by reverse sorting the input array. If there are still more than IDBLIM
00135        above the threshold, then revise the thresold until there aren't. Then
00136        use the variable npl2 to stop the rest of the routine accessing any of
00137        the fainter pixels in the list. This is to restrict processing time
00138        for large extended objects */
00139  
00140     curthr = smul*oldthr;
00141     sort_on_zsm_rev(npl,pl);
00142     while (1) {
00143         npl2 = 0;
00144         while (pl[npl2].zsm > curthr && npl2 < npl-1)
00145             npl2++;    
00146         if (npl2 > IDBLIM) 
00147             curthr += oldthr;
00148         else
00149             break;
00150     }
00151 
00152     /* If there are fewer pixels above the new threshold than the minimum 
00153        specified in the input parameters, then you have no reason to be here */
00154 
00155     if (npl2 < ipix) {
00156         *nbit = 1;
00157         return;
00158     }
00159 
00160     /* Get a new ap structure */
00161 
00162     ap2.lsiz = ap->lsiz;
00163     ap2.csiz = ap->csiz;
00164     ap2.multiply = 1;
00165     ap2.ipnop = ipixo2;
00166     ap2.areal_offset = offset;
00167     ap2.fconst = fconst;
00168     ap2.mflag = cpl_calloc((ap2.lsiz)*(ap2.csiz),sizeof(unsigned char*));
00169 
00170     /* Main analysis loop at new thresholds */
00171 
00172     *nbit = 0;
00173     nbitprev = 0;
00174     apinit(&ap2);
00175     while (1) {
00176         nexthr = MAX(curthr+oldthr,curthr*tmul);
00177         
00178         /* Locate objects in this cluster */
00179 
00180         ap2.thresh = curthr;
00181         apclust(&ap2,npl2,pl);
00182         check_term(&ap2,&nobj,results,ipks,&toomany);
00183         apreinit(&ap2);
00184         if (nobj == 0)
00185             break;
00186 
00187         /* For each image check the number of points above the next threshold
00188            and flag. Do a moments analysis of each object */
00189 
00190         for (i = 0; i < nobj; i++) {
00191 
00192             /* Ok, has this object already been detected?  If so, then
00193                load the new results into the parmnew array. We'll check
00194                whether it's worthwhile updating the master parameters
00195                list later */
00196 
00197             isnew = 1;
00198             xb = results[i][1];
00199             yb = results[i][2];
00200             sxx = MAX(1.0,results[i][4]);
00201             syy = MAX(1.0,results[i][6]);
00202             for (k = 0; k < nbitprev; k++) {
00203                 dx = xb - parm[k][1];
00204                 dy = yb - parm[k][2];
00205                 radius2 = dx*dx/sxx + dy*dy/syy;
00206                 if ((ibitx[k] == ipks[i][0] && ibity[k] == ipks[i][1]) || 
00207                     radius2 < 1.0) {
00208                     isnew = 0;
00209                     for (kk = 0; kk < NPAR; kk++)
00210                         parmnew[k][kk] = results[i][kk];
00211                     ibitl[k] = (int)results[i][NPAR];
00212                     break;
00213                 }
00214             }
00215 
00216             /* If this is a new one and it's above the minimum threshold 
00217                then check to make sure you don't already have too many. 
00218                If you do, then flag this and break out to the next iteration.
00219                If there is room for another, then store the moments analysis
00220                profile */
00221 
00222             if (isnew && results[i][0] > xintmn) {
00223                 if (*nbit >= IMNUM) {
00224                     *nbit = IMNUM;
00225                     toomany = 1;
00226                     break;
00227                 }
00228                 ibitx[*nbit] = ipks[i][0];
00229                 ibity[*nbit] = ipks[i][1];
00230                 for (kk = 0; kk < NPAR; kk++)
00231                     parm[*nbit][kk] = results[i][kk];
00232                 ibitl[*nbit] = (int)results[i][NPAR];
00233                 (*nbit)++;
00234             }
00235         } /* End of object loop */
00236 
00237         /* If too many objects were found, then skip the next bit...otherwise
00238            go through and update parameters if necessary. This block of
00239            code is a bit of a place holder waiting for something better to
00240            be worked out...*/
00241 
00242         if (! toomany) {
00243             if (*nbit > nbitprev && nbitprev > 0) {
00244                 for (i = 0; i < nbitprev; i++) 
00245                     iupdate[i] = 0;
00246                 for (j = nbitprev; j < *nbit; j++) {
00247                     distmax = 0.0;
00248                     iwas = 0;
00249                     for (i = 0; i < nbitprev; i++) {
00250                         if (parmnew[i][0] > 0.0) {
00251                             dx = parmnew[i][1] - parm[i][1];
00252                             dy = parmnew[i][2] - parm[i][2];
00253                             radius2 = dx*dx + dy*dy;
00254                             if (radius2 > distmax) {
00255                                 iwas = i;
00256                                 distmax = radius2;
00257                             }
00258                         }
00259                     }
00260                     iupdate[iwas] = 1;
00261                 }
00262                 for (i = 0; i < nbitprev; i++)
00263                     if (iupdate[i] == 1 && parmnew[i][0] > 0.0) 
00264                         for (j = 0; j < NPAR; j++) 
00265                             parm[i][j] = parmnew[i][j];
00266             }
00267 
00268             /* Reset the update flag and prepare for next iteration*/
00269 
00270             for (i = 0; i <= *nbit; i++)
00271                 parmnew[i][0] = -1.0;
00272             nbitprev = *nbit;
00273         }
00274 
00275         /* Where do we cut in the list now? */
00276 
00277 
00278         npl3 = 0;
00279         while (pl[npl3].zsm > nexthr && npl3 < npl2-1)
00280             npl3++;    
00281         npl2 = npl3;
00282 
00283         /* Now, do we need to move onto the next threshold? */
00284 
00285         if (npl2 == 0 || toomany || nexthr >= itmaxlim) 
00286             break;
00287 
00288         /* If so, then reset some variables and continue */
00289 
00290         curthr = nexthr;
00291 
00292     } /* End of main analysis loop */
00293 
00294     /* Free up some workspace */
00295 
00296     free(ap2.mflag);
00297     apclose(&ap2);
00298 
00299     /* If there is only one then we can exit now */
00300 
00301     if (*nbit == 1)
00302         return;
00303 
00304     /* Find out which images terminated properly and remove those that didn't */
00305     j = -1;
00306     for (k = 0; k < *nbit; k++) {
00307         /* Commented this out as checking for terminations seems to miss some
00308            if the total flux above the threshold for an object is negative */
00309 
00310 /*         if (ibitl[k] == 1 && parm[k][0] > xintmn) { */ 
00311         if (parm[k][0] > xintmn) {
00312             j++;
00313             if (j != k)
00314                 for (i = 0; i < NPAR; i++)
00315                     parm[j][i] = parm[k][i];
00316         }
00317     }
00318     *nbit = j + 1;
00319     for (j = 0; j < *nbit; j++) {
00320         bitx[j] = 0.0;
00321         bitl[j] = 0.0;
00322     }
00323 
00324     /* For each image find true areal profile levels and iterate to find 
00325        local continuum */
00326 
00327     iter = 0;
00328     sumint = 0.0;
00329     lastone = 0;
00330     while (iter < NITER) {
00331         iter++;
00332         
00333         /* Loop for each of the objects and create a level vs radius array */
00334 
00335         for (k = 0; k < *nbit; k++) {
00336             if (parm[k][0] < 0.0)
00337                 continue;
00338             xlevol = logf(parm[k][7] + parm[k][3] - bitl[k]); /* Pk + newthresh - cont */
00339             xlevel = xlevol;
00340             radold = 0.0;
00341             radius = radold;
00342             slope = 1.0;
00343             ic = 0;
00344             for (i = 1; i <= NAREAL; i++) {
00345                 jj = NPAR - i;
00346                 ii = NAREAL - i;
00347                 xx = (float)ii + offset;
00348                 if (parm[k][jj] > 0.5) {
00349                     if (ii == 0) 
00350                         xlevel = logf(parm[k][3] - bitl[k] + 0.5);
00351                     else 
00352                         xlevel = logf(powf(2.0,xx) - oldthr + parm[k][3] -
00353                                       bitl[k] - 0.5);
00354                     radius = sqrt(parm[k][jj]/CPL_MATH_PI);
00355                     xdat[ic] = xlevel;
00356                     xcor[ic++] = radius;
00357                     dlbydr = (xlevol - xlevel)/MAX(0.01,radius-radold);
00358                     wt = MIN(1.0,MAX((radius-radold)*5.0,0.1));
00359                     slope = (1.0 - 0.5*wt)*slope + 0.5*wt*MIN(5.0,dlbydr);
00360                     radold = radius;
00361                     xlevol = xlevel;
00362                 }
00363             }
00364 
00365             /* If this is not the last iteration then work out the effect
00366                on the local continuum from each image */
00367 
00368             if (! lastone) {
00369                 for (i = 0; i < *nbit; i++) {
00370                     if (parm[i][0] >= 0.0 && i != k) {
00371                         dx = parm[k][1] - parm[i][1];
00372                         dy = parm[k][2] - parm[i][2];
00373                         dist = sqrtf(dx*dx + dy*dy);
00374                         xeff = xlevel - MAX(0.0,MIN(50.0,slope*(dist-radius)));
00375                         bitx[i] += expf(xeff);
00376                     }
00377                 }
00378 
00379             /* If this is the last iteration loop, then update the parameters
00380                before exiting*/
00381             
00382             } else {
00383                 if (ic > 2) {
00384                     polynm(xdat,xcor,ic,polycf,3,0);
00385                     ttt = polycf[1] + 2.0*polycf[2]*radius;
00386                 } else
00387                     ttt = 0.0;
00388                 slope = MAX(0.1,MAX(-ttt,slope));
00389                 radthr = radius + (xlevel - algthr)/slope;
00390                 if (radthr > radmax) {
00391                     slope = 1.0;
00392                     radthr = radmax;
00393                 }
00394                 
00395                 /* Pixel area */
00396                
00397                 delb = parm[k][8]*(parm[k][3] - bitl[k]);
00398                 parm[k][8] = CPL_MATH_PI*radthr*radthr;
00399 
00400                 /* Peak height */
00401 
00402                 parm[k][7] += (parm[k][3] - bitl[k]);
00403                 
00404                 /* Intensity */
00405 
00406                 deli = 2.0*CPL_MATH_PI*((parm[k][3] - bitl[k])*(1.0 + slope*radius) -
00407                                  oldthr*(1.0 + slope*radthr))/(slope*slope);
00408                 parm[k][0] += delb + MAX(0.0,deli);
00409                 for (i = 0; i < 7; i++)
00410                     parm[k][i+9] = -1.0;
00411                 if (parm[k][0] > xintmn)
00412                     sumint += parm[k][0];
00413             }
00414         }
00415 
00416         /* If this is not the last iteration then check and see how the
00417            continuum estimates are converging. If they appear to be converging
00418            then let the next iteration be the last one. */
00419 
00420         if (! lastone) {
00421             conv = 1;
00422             for (i = 0; i < *nbit; i++) {
00423                 if (parm[i][0] >= 0.0) {
00424                     if (abs(bitx[i] - bitl[i]) > 3.0)
00425                         conv = 0;
00426                     bitl[i] = bitx[i];
00427                     bitx[i] = 0;
00428                     bitl[i] = MIN(bitl[i],NINT(parm[i][3]-oldthr));
00429                 }
00430             }
00431             lastone = (conv || (iter == NITER-1));
00432         } else {
00433             break;
00434         }
00435     }
00436 
00437     /* Find the scaling if needed */
00438 
00439     if (sumint == 0.0) {
00440         *nbit = 1;
00441         return;
00442     } else 
00443         ratio = total/sumint;
00444     for (i = 0; i < *nbit; i++)
00445         parm[i][0] = ratio*parm[i][0];
00446 }
00447 
00448 /*---------------------------------------------------------------------------*/
00468 /*---------------------------------------------------------------------------*/
00469         
00470 static void sort_on_zsm_rev(int npts, plstruct *pts) {
00471     int i,j,ii,jj,ifin;
00472     plstruct tmp;
00473 
00474     jj = 4;
00475     while (jj < npts)
00476         jj = 2*jj;
00477     jj = MIN(npts,(3*jj)/4-1);
00478     while (jj > 1) {
00479         jj = jj/2;
00480         ifin = npts - jj;
00481         for (ii = 0; ii < ifin; ii++) {
00482             i = ii;
00483             j = i + jj;
00484             if (pts[i].zsm > pts[j].zsm)
00485                 continue;
00486             tmp = pts[j];
00487             do {
00488                 pts[j] = pts[i];
00489                 j = i;
00490                 i = i - jj;
00491                 if (i < 0) 
00492                     break;
00493             } while (pts[i].zsm <= tmp.zsm);
00494             pts[j] = tmp;
00495         }
00496     }
00497 }
00498 
00499 /*---------------------------------------------------------------------------*/
00521 /*---------------------------------------------------------------------------*/
00522 
00523 static void moments_thr(ap_t *ap, float results[NPAR+1], int ipk[2]) {
00524     int i,np,nnext;
00525     float x,y,xoff,yoff,xsum,ysum,xsumsq,ysumsq,tsum,xysum,t,tmax,twelfth;
00526     float xbar,ybar,sxx,syy,sxy,fconst,offset,xsum_w,ysum_w,wsum,w;
00527     plstruct *plarray;
00528 
00529     /* Copy some stuff to local variables */
00530 
00531     fconst = ap->fconst;
00532     offset = ap->areal_offset;
00533     plarray = ap->plarray;
00534     np = ap->npl_pix;
00535 
00536     /* Initialise a few things */
00537 
00538     xoff = xbar_start;
00539     yoff = ybar_start;
00540     xsum = 0.0;
00541     ysum = 0.0;
00542     xsum_w = 0.0;
00543     ysum_w = 0.0;
00544     wsum = 0.0;
00545     xsumsq = 0.0;
00546     ysumsq = 0.0;
00547     tsum = 0.0;
00548     xysum = 0.0;
00549     tmax = plarray[0].z - curthr;
00550     ipk[0] = plarray[0].x;
00551     ipk[1] = plarray[0].y;
00552     twelfth = 1.0/12.0;
00553     for (i = 8; i < NPAR; i++)
00554         results[i] = 0.0;
00555 
00556     /* Do a moments analysis on an object */
00557 
00558     nnext = 0;
00559     for (i = 0; i < np; i++) {
00560         x = (float)plarray[i].x - xoff;
00561         y = (float)plarray[i].y - yoff;
00562         t = plarray[i].z - curthr;
00563         w = plarray[i].zsm - curthr;
00564         if (w > nexthr)
00565             nnext++;
00566         xsum += t*x;
00567         ysum += t*y;
00568         tsum += t;
00569         xsum_w += w*t*x;
00570         ysum_w += w*t*y;
00571         wsum += w*t;
00572         xsumsq += (x*x + twelfth)*t;
00573         ysumsq += (y*y + twelfth)*t;
00574         xysum += x*y*t;
00575         update_ov(results+8,t,oldthr,fconst,offset);
00576         if (t > tmax) {
00577             ipk[0] = plarray[i].x;
00578             ipk[1] = plarray[i].y;
00579             tmax = t;
00580         }
00581     }
00582 
00583     /* Check that the total intensity is enough and if it is, then do
00584        the final results. Use negative total counts to signal an error */
00585 
00586     if (tsum > 0.0) {
00587         results[0] = tsum;
00588     } else {
00589         results[0] = -1.0;
00590         tsum = 1.0;
00591     }
00592     xbar = xsum/tsum;
00593     ybar = ysum/tsum;
00594     sxx = MAX(0.0,(xsumsq/tsum-xbar*xbar));
00595     syy = MAX(0.0,(ysumsq/tsum-ybar*ybar));
00596     sxy = xysum/tsum - xbar*ybar;
00597     wsum = MAX(1.0,wsum);
00598     xbar = xsum_w/wsum;
00599     ybar = ysum_w/wsum;
00600     xbar += xoff;
00601     ybar += yoff;
00602     xbar = MAX(1.0,MIN(xbar,ap->lsiz));
00603     ybar = MAX(1.0,MIN(ybar,ap->csiz));
00604 
00605     /* Store the results now */
00606 
00607     results[1] = xbar;
00608     results[2] = ybar;
00609     results[3] = curthr;
00610     results[4] = sxx;
00611     results[5] = sxy;
00612     results[6] = syy;
00613     results[7] = tmax;
00614     results[NPAR] = ((nnext > ap->ipnop && nexthr < lasthr) ? 0 : 1);
00615 }
00616 
00617 /*---------------------------------------------------------------------------*/
00642 /*---------------------------------------------------------------------------*/
00643         
00644 static void update_ov(float iap[NAREAL], float t, float thresh, float fconst, 
00645                       float offset) {
00646     int nup,i;
00647 
00648     /* Get out of here if the intensity is too small */
00649 
00650     if (t <= 0.0) 
00651         return;
00652 
00653     /* Otherwise update the relevant profile counts */
00654 
00655     nup = MAX(1,MIN(NAREAL,(int)(logf(t+thresh)*fconst-offset)+1));
00656     for (i = 0; i < nup; i++)
00657         iap[i] += 1.0;
00658 }
00659 
00660 /*---------------------------------------------------------------------------*/
00686 /*---------------------------------------------------------------------------*/
00687 
00688 static void check_term(ap_t *ap, int *nobj, float parm[IMNUM][NPAR+1], 
00689                        int peaks[IMNUM][2], int *toomany) {
00690     int ip,i,ipks[2];
00691     float momresults[NPAR+1];
00692 
00693     /* Search through all possible parents */
00694 
00695     *nobj = 0;
00696     *toomany = 0;
00697     for (ip = 1; ip <= ap->maxip; ip++) {
00698         if (ap->parent[ip].pnop != -1) {
00699 /*             if (ap->parent[ip].pnop == ap->parent[ip].growing) { */
00700   
00701                 /* That's a termination: */
00702  
00703                 if ((ap->parent[ip].pnop >= ap->ipnop &&
00704                      ap->parent[ip].touch == 0)) {
00705                     extract_data(ap,ip);
00706                     moments_thr(ap,momresults,ipks);
00707                     if (momresults[0] > 0.0) {
00708                         if (*nobj == IMNUM-1) {
00709                             *toomany = 1;
00710                             break;
00711                         }
00712                         for (i = 0; i <= NPAR; i++) 
00713                             parm[*nobj][i] = momresults[i];
00714                         for (i = 0; i < 2; i++)
00715                             peaks[*nobj][i] = ipks[i];
00716                         (*nobj)++;
00717                     }
00718                 }
00719                 restack(ap,ip);
00720 /*          } else { */
00721 
00722 /*              /\* This parent still active: *\/ */
00723  
00724 /*              ap->parent[ip].growing = ap->parent[ip].pnop; */
00725 /*          } */
00726         }
00727     }
00728 }
00729 
00732 /*
00733 
00734 $Log: imcore_overlp.c,v $
00735 Revision 1.9  2010/10/05 09:11:33  jim
00736 Modified moments_thr to remove possibility of uninitialised values
00737 
00738 Revision 1.8  2010/09/09 12:09:57  jim
00739 Added docs
00740 
00741 Revision 1.7  2009/09/09 09:42:25  jim
00742 modified to use CPL defined macros for constants
00743 
00744 Revision 1.6  2009/01/23 12:24:33  jim
00745 Fixed bugs in pixel flagging
00746 
00747 Revision 1.5  2007/12/19 13:17:16  jim
00748 Fixed small bug affecting how parameters get updated
00749 
00750 Revision 1.4  2006/08/21 09:07:29  jim
00751 Modified centring algorithm
00752 
00753 Revision 1.3  2006/06/29 13:43:05  jim
00754 Modified to take some stuff out of integer arithmetic
00755 
00756 Revision 1.2  2006/04/20 11:15:19  jim
00757 A few minor bugs fixed
00758 
00759 Revision 1.1  2005/09/13 13:25:29  jim
00760 Initial entry after modifications to make cpl compliant
00761 
00762 
00763 */

Generated on 7 Feb 2011 for VIRCAM Pipeline by  doxygen 1.6.1