imcore_phopt.c

00001 /* $Id: imcore_phopt.c,v 1.6 2010/09/09 12:09:57 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/09/09 12:09:57 $
00024  * $Revision: 1.6 $
00025  * $Name: vcam-1_3_0 $
00026  */
00027 
00028 #include <stdio.h>
00029 #include <stdlib.h>
00030 #include <math.h>
00031 
00032 #include <cpl.h>
00033 
00034 #include "floatmath.h"
00035 #include "util.h"
00036 #include "imcore.h"
00037 #include "ap.h"
00038 
00039 /* Function Prototypes */
00040 
00041 static void dchole (double a[IMNUM+1][IMNUM+1], double b[IMNUM+1], int n);
00042 static float fraction (float x, float y, float r_out);
00043 
00044 /* does multiple profile fitting to determine intensities */
00045 
00048 /*---------------------------------------------------------------------------*/
00089 /*---------------------------------------------------------------------------*/
00090 
00091 extern void phopt(ap_t *ap, float parm[IMNUM][NPAR], int nbit, int naper, 
00092                   float apertures[], float cflux[], float badpix[], 
00093                   int nrcore) {
00094     double aa[IMNUM+1][IMNUM+1],bb[IMNUM+1];
00095     float d,arg,*map,rcirc;
00096     float cn,parrad,xmin,xmax,xi,yi,ymin,ymax;
00097     float t,xj,yj,cnsq,tj,xk,yk,tk;
00098     int i,ii,j,kk,ix1,ix2,iy1,iy2,nx,ny,k,iaper;
00099     unsigned char *mflag,mf;
00100 
00101     /* Set up some local variables */
00102 
00103     map = ap->indata;
00104     mflag = ap->mflag;
00105     nx = ap->lsiz;
00106     ny = ap->csiz;   
00107 
00108     /* Loop for each of the apertures */
00109 
00110     for (iaper = 0; iaper < naper; iaper++) {
00111         rcirc = apertures[iaper];
00112         parrad = rcirc + 0.5;
00113         cn = 1.0/(CPL_MATH_PI*rcirc*rcirc);  /* profile normalising constant */
00114         cnsq = cn*cn;
00115     
00116         /* set up covariance matrix - analytic special case for cores */
00117 
00118         for(i = 0; i < nbit; i++) {
00119             aa[i][i] = cn;                 /* overlaps totally area=pi*r**2 */
00120             if(nbit > 1) {
00121                 xi = parm[i][1];
00122                 yi = parm[i][2];
00123                 for(j = i+1; j < nbit; j++) {
00124                     d = sqrtf((xi-parm[j][1])*(xi-parm[j][1]) 
00125                         + (yi-parm[j][2])*(yi-parm[j][2]));
00126                     if(d >= 2.0*rcirc) {
00127                         aa[j][i] = 0.0;
00128                         aa[i][j] = aa[j][i];
00129                     } else {
00130                         arg = d/(2.0*rcirc);
00131                         aa[j][i] = cnsq*2.0*rcirc*rcirc*
00132                             (acosf(arg)-arg*(sqrtf(1.0-arg*arg)));
00133                         aa[i][j] = aa[j][i];
00134                     }
00135                 }
00136             }
00137         }
00138 
00139         /* clear accumulators */
00140 
00141         for(i = 0; i < nbit; i++) 
00142             bb[i] = 0.0;
00143 
00144         /* generate image-blend outer boundaries */
00145 
00146         xmin = 1.0e6;
00147         xmax = -1.0e6;
00148         ymin = 1.0e6;
00149         ymax = -1.0e6;
00150         for(i = 0; i < nbit; i++) {
00151             xi = parm[i][1];
00152             yi = parm[i][2];
00153             xmin = MIN(xmin, xi);
00154             xmax = MAX(xmax, xi);
00155             ymin = MIN(ymin, yi);
00156             ymax = MAX(ymax, yi);
00157         }
00158         ix1 = MAX(0,(int)(xmin-parrad)-1);
00159         ix2 = MIN(nx-1,(int)(xmax+parrad));
00160         iy1 = MAX(0,(int)(ymin-parrad)-1);
00161         iy2 = MIN(ny-1,(int)(ymax+parrad));
00162 
00163         /* now go through pixel region */
00164 
00165         for(ii = iy1; ii <= iy2; ii++) {
00166             kk = ii*nx;
00167             for(i = ix1; i <= ix2; i++) {
00168                 mf = mflag[kk+i];
00169                 if (mf == MF_ZEROCONF || mf == MF_STUPID_VALUE) {
00170                     for (j = 0; j < nbit; j++) {
00171                         xj = i - parm[j][1] + 1.0;
00172                         yj = ii - parm[j][2] + 1.0;
00173                         tj = fraction(xj,yj,rcirc);
00174                         aa[j][j] -= tj*tj*cnsq;
00175                         for (k = j + 1; k < nbit; k++) {
00176                             xk = i - parm[k][1] + 1.0;
00177                             yk = ii - parm[k][2] + 1.0;
00178                             tk = fraction(xk,yk,rcirc);
00179                             aa[k][j] -= tk*tj*cnsq;
00180                             aa[j][k] = aa[k][j];
00181                         }
00182                         if (iaper == nrcore)
00183                             badpix[j] += tj;
00184                     }
00185                 } else if (mf == MF_CLEANPIX || mf == MF_OBJPIX ||
00186                            mf == MF_SATURATED) {
00187                     t = map[kk+i];   
00188                     for(j = 0; j < nbit; j++) {
00189                         xj = i - parm[j][1] + 1.0;
00190                         yj = ii - parm[j][2] + 1.0;
00191                         bb[j] += fraction(xj,yj,rcirc)*t;
00192                     }
00193                 } 
00194             }
00195         }
00196 
00197         /* Trivial solution for single object */
00198 
00199         if (nbit == 1) {
00200             cflux[iaper] = bb[0];
00201 
00202         /* solve for profile intensities */
00203 
00204         } else {
00205             for (i = 0; i < nbit; i++) 
00206                 aa[i][i] = MAX(aa[i][i],cnsq);
00207             dchole(aa,bb,nbit);
00208             for(i = 0; i < nbit; i++) 
00209                 cflux[i*naper+iaper] = cn*bb[i];
00210         }
00211     }
00212 }
00213 
00216 /* CHOLEsky decomposition of +ve definite symmetric matrix to solve Ax = b */
00217 
00218 static void dchole (double a[IMNUM+1][IMNUM+1], double b[IMNUM+1], int n) {
00219   double sum, l[IMNUM+1][IMNUM+1], y[IMNUM+1];
00220   double aveigv, offset;
00221   int i, j, k;
00222 
00223 restart:
00224     l[0][0] = sqrt(a[0][0]);
00225 
00226     for(k = 1; k < n; k++) {
00227         for(j = 0; j <= k-1; j++) {
00228             sum = a[j][k];
00229             if(j != 0) 
00230                 for(i = 0; i <= j-1; i++) 
00231                     sum -= l[i][k]*l[i][j];
00232             l[j][k] = sum/l[j][j];
00233         }
00234         sum = a[k][k];
00235         for(i = 0; i <= k-1; i++) 
00236             sum -= l[i][k]*l[i][k];
00237         if(sum <= 0.0) {
00238 /*          fprintf(stderr, "dchole: warning: matrix ill-conditioned\n"); */
00239             aveigv = a[0][0];
00240             for(i = 1; i < n; i++) 
00241                 aveigv += a[i][i];
00242             /* max eigenvalue < trace */
00243             offset = 0.1*aveigv/((double) n);
00244             for(i = 0; i < n; i++) 
00245                 a[i][i] += offset;
00246 /*          fprintf(stderr, "dchole: Offset added to diagonal = %f\n", offset); */
00247             goto restart;
00248         }
00249         l[k][k] = sqrt(sum);
00250     }
00251 
00252     /* solve Ly = b */
00253 
00254     y[0] = b[0]/l[0][0];
00255     for(i = 1; i < n; i++) {
00256         sum = b[i];
00257         for(k = 0; k <= i-1; k++) 
00258             sum -= l[k][i]*y[k];
00259         y[i] = sum/l[i][i];
00260     }
00261 
00262     /* solve L(T)x = y */
00263 
00264     b[n-1] = y[n-1]/l[n-1][n-1];
00265     for(i = n-2; i >= 0; i--) {
00266         sum = y[i];
00267         for(k = i+1; k < n; k++) 
00268             sum -= l[i][k]*b[k];
00269         b[i] = sum/l[i][i];
00270     }
00271 }
00272 
00273 /* returns fraction of pixel bounded by 0 -  r_out
00274  * x,y coordinates relative to centre
00275  * Uses linear approximation ok if pixel located >>1 away from centre */
00276 
00277 static float fraction (float x, float y, float r_out) {
00278     float r,t,x_a,x_b,frac,tanao2,cosa,tanp2a,sqrt2o2;
00279 
00280     r = sqrtf(x*x + y*y);
00281     sqrt2o2 = 0.5*CPL_MATH_SQRT2;
00282 
00283     /* is it worth bothering? */
00284 
00285     if(r > r_out+sqrt2o2) 
00286         return(0.0);
00287 
00288     /* is it trivially all in? */
00289 
00290     if(r < r_out-sqrt2o2) 
00291         return(1.0);
00292 
00293     /* bugger - have to do some work then ... ok first ...
00294      * use 8-fold symmetry to convert to 0-45 degree range */
00295 
00296     x = fabsf(x);
00297     y = fabsf(y);
00298     if(y > x) {
00299         t = x;
00300         x = y;
00301         y = t;
00302     }
00303 
00304     /* If the angles are too close to cardinal points, then fudge something */
00305 
00306     if (x > 0.0 && y > 0.0) {
00307         tanao2 = 0.5*y/x;
00308         tanp2a = x/y;
00309         cosa = x/sqrt(x*x + y*y);
00310     } else {
00311         tanao2 = 0.00005;
00312         tanp2a = 10000.0;
00313         cosa = 1.0;
00314     }
00315 
00316     /* only outer radius - compute linear intersections top and bot of pixel */
00317 
00318     x_a = x - tanao2 + (r_out - r)/cosa;
00319     if(x_a < x+0.5) {
00320 
00321         /* intersects */
00322 
00323         x_b = x + tanao2 + (r_out - r)/cosa;
00324 
00325         /* three cases to consider */
00326 
00327         if(x_a < x-0.5)
00328             frac = 0.5*MAX(0.0,x_b-(x-0.5))*MAX(0.0,x_b-(x-0.5))*tanp2a;
00329         else {
00330             if(x_b > x+0.5)
00331                 frac = 1.0 - 0.5*(x+0.5-x_a)*(x+0.5-x_a)*tanp2a;
00332             else
00333                 frac = 0.5-(x-x_a)+0.5*(x_b-x_a);
00334         }
00335     } else  /* missed entirely */
00336         frac = 1.0;
00337 
00338     return(frac);
00339 }
00340 
00341 /*
00342 
00343 $Log: imcore_phopt.c,v $
00344 Revision 1.6  2010/09/09 12:09:57  jim
00345 Added docs
00346 
00347 Revision 1.5  2009/09/09 09:42:25  jim
00348 modified to use CPL defined macros for constants
00349 
00350 Revision 1.4  2009/01/23 12:24:33  jim
00351 Fixed bugs in pixel flagging
00352 
00353 Revision 1.3  2009/01/20 09:31:06  jim
00354 Fixed problem with pixel flagging
00355 
00356 Revision 1.2  2008/10/13 08:11:19  jim
00357 Fixed pixel masking scheme
00358 
00359 Revision 1.1  2005/09/13 13:25:29  jim
00360 Initial entry after modifications to make cpl compliant
00361 
00362 
00363 */

Generated on 15 Mar 2012 for VIRCAM Pipeline by  doxygen 1.6.1