VIRCAM Pipeline  1.3.3
vircam/vircam_imstack.c
1 /* $Id: vircam_imstack.c,v 1.16 2013-10-15 16:46:13 jim Exp $
2  *
3  * This file is part of the VIRCAM Pipeline
4  * Copyright (C) 2005 Cambridge Astronomy Survey Unit
5  *
6  * This program is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2 of the License, or
9  * (at your option) any later version.
10  *
11  * This program is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License
17  * along with this program; if not, write to the Free Software
18  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19  */
20 
21 /*
22  * $Author: jim $
23  * $Date: 2013-10-15 16:46:13 $
24  * $Revision: 1.16 $
25  * $Name: not supported by cvs2svn $
26  */
27 
28 /* Includes */
29 
30 #ifdef HAVE_CONFIG_H
31 #include <config.h>
32 #endif
33 
34 #include <math.h>
35 #include <cpl.h>
36 #include "vircam_mods.h"
37 #include "vircam_utils.h"
38 #include "vircam_fits.h"
39 #include "vircam_stats.h"
40 #include "vircam_pfits.h"
41 #include "vircam_wcsutils.h"
42 #include "vircam_tfits.h"
43 
44 typedef struct {
45  double crval[2];
46  double crpix[2];
47  double cd[4];
48  double secd;
49  double tand;
50  double pv21;
51  double pv23;
52  double pv25;
53  long nx;
54  long ny;
55 } mywcs;
56 
57 typedef struct {
58  int nx;
59  int ny;
60  vir_fits *fname;
61  vir_fits *conf;
62  float *data;
63  int *cdata;
64  cpl_binary *bpm;
65  mywcs *vwcs;
66  float xoff;
67  float yoff;
68  cpl_array *trans;
69  float sky;
70  float skydiff;
71  float noise;
72  float expscale;
73  float weight;
74 } dstrct;
75 
76 static dstrct *fileptrs = NULL;
77 static vir_fits *outfits = NULL;
78 static vir_fits *outconf = NULL;
79 static mywcs *outwcs = NULL;
80 
81 static float lsig;
82 static float hsig;
83 static int nxo;
84 static int nyo;
85 static float *outdata;
86 static int *outdatac;
87 static int nim;
88 static float sumweight;
89 
90 static void outloc(mywcs *win, double xin, double yin, mywcs *wout,
91  cpl_array *trans, double *xout, double *yout);
92 static void seeing_wt(void);
93 static void backgrnd_ov(void);
94 static void skyest(float *data, int *cdata, long npts, float thresh,
95  float *skymed, float *skynoise);
96 static cpl_array *transinit(void);
97 static void fileptrs_close(void);
98 static cpl_table *mknewtab(vir_tfits *catref, dstrct *dd);
99 static void output_files(void);
100 static int stack_nn(void);
101 static int stack_lin(void);
102 static void do_averages(int ncontrib, float *data, float *wconf, float *conf,
103  unsigned char *id, float lclip, float hclip,
104  float extra, short int *lowcl, short int *highcl,
105  float *outval, float *outvalc);
106 static mywcs *getwcsinfo(cpl_propertylist *plist);
107 static void diffxy(mywcs *wref, mywcs *wprog, float *dx, float *dy);
108 static void tidy(void);
109 
112 /*---------------------------------------------------------------------------*/
172 /*---------------------------------------------------------------------------*/
173 
174 extern int vircam_imstack(vir_fits **inf, vir_fits **inconf, vir_tfits **cats,
175  int nimages, int nconfs, float lthr, float hthr,
176  int method, int seeing, vir_fits **out,
177  vir_fits **outc, int *status) {
178 
179  int i,nm;
180  dstrct *dd;
181  float expref,exposure,xoff,yoff;
182  cpl_matrix *in;
183  vir_tfits *catref;
184  cpl_table *newtab,*mtab;
185  const char *fctid = "vircam_imstack";
186 
187  /* Inherited status */
188 
189  *out = NULL;
190  *outc = NULL;
191  if (*status != VIR_OK)
192  return(*status);
193 
194  /* Is there any point in being here? */
195 
196  if (nimages == 0) {
197  cpl_msg_error(fctid,"No input files to combine");
198  tidy();
199  FATAL_ERROR
200  }
201 
202  /* Initialise some global variables */
203 
204  lsig = lthr;
205  hsig = hthr;
206  nim = nimages;
207 
208  /* Initialise a useful matrix */
209 
210  in = cpl_matrix_new(1,2);
211  cpl_matrix_set(in,0,0,0.0);
212  cpl_matrix_set(in,0,1,0.0);
213 
214  /* Allocate file struct array and fill in some values */
215 
216  fileptrs = cpl_malloc(nimages*sizeof(dstrct));
217  (void)vircam_pfits_get_exptime(vircam_fits_get_phu(inf[0]),&exposure);
218  expref = max(0.5,exposure);
219  catref = NULL;
220  for (i = 0; i < nimages; i++) {
221  dd = fileptrs + i;
222  dd->fname = inf[i];
223  dd->data = cpl_image_get_data_float(vircam_fits_get_image(inf[i]));
224  if (nconfs == 0) {
225  dd->conf = NULL;
226  } else if (nconfs == 1) {
227  dd->conf = inconf[0];
228  dd->cdata = cpl_image_get_data_int(vircam_fits_get_image(inconf[0]));
229  } else {
230  dd->conf = inconf[i];
231  dd->cdata = cpl_image_get_data_int(vircam_fits_get_image(inconf[i]));
232  }
233  (void)vircam_pfits_get_exptime(vircam_fits_get_phu(dd->fname),&exposure);
234  dd->expscale = exposure/expref;
235 
236  /* Zero the bad pixel mask */
237 
238  cpl_image_accept_all(vircam_fits_get_image(inf[i]));
239  dd->bpm = cpl_mask_get_data(cpl_image_get_bpm(vircam_fits_get_image(inf[i])));
240 
241  /* Read the WCS from the header */
242 
243  dd->vwcs = getwcsinfo(vircam_fits_get_ehu(inf[i]));
244 
245  /* Get the image size */
246 
247  dd->nx = (int)cpl_image_get_size_x(vircam_fits_get_image(dd->fname));
248  dd->ny = (int)cpl_image_get_size_y(vircam_fits_get_image(dd->fname));
249 
250  /* Double check to make sure the confidence maps and images have the
251  same dimensions */
252 
253  if ((int)cpl_image_get_size_x(vircam_fits_get_image(dd->conf)) != dd->nx ||
254  (int)cpl_image_get_size_y(vircam_fits_get_image(dd->conf)) != dd->ny) {
255  cpl_msg_error(fctid,"Image %s and Confidence map %s don't match",
256  vircam_fits_get_fullname(dd->fname),
257  vircam_fits_get_fullname(dd->conf));
258  cpl_matrix_delete(in);
259  tidy();
260  FATAL_ERROR
261  }
262 
263  /* Get a rough shift between the frames and the first one */
264 
265  diffxy(fileptrs->vwcs,dd->vwcs,&(dd->xoff),&(dd->yoff));
266 
267  /* Do we have catalogues? If so, then work out the transformation
268  corrections from them */
269 
270  if (cats != NULL) {
271  if (catref == NULL && cats[i] != NULL)
272  catref = cats[i];
273  if (i == 0 || cats[i] == NULL) {
274  dd->trans = transinit();
275  } else {
276  newtab = mknewtab(catref,dd);
277  (void)vircam_matchxy(vircam_tfits_get_table(cats[i]),newtab,
278  1024,&xoff,&yoff,&nm,&mtab,status);
279  if (*status != VIR_OK) {
280  freearray(dd->trans);
281  dd->trans = transinit();
282  *status = VIR_OK;
283  } else {
284  (void)vircam_platexy(mtab,6,&(dd->trans),status);
285  }
286  freetable(mtab);
287  freetable(newtab);
288  }
289  } else {
290  dd->trans = NULL;
291  }
292  }
293  cpl_matrix_delete(in);
294 
295  /* Get the background levels in the overlap regions */
296 
297  backgrnd_ov();
298 
299  /* Do the seeing weighting if you want to */
300 
301  if (seeing)
302  seeing_wt();
303 
304  /* Create the output file */
305 
306  output_files();
307 
308  /* Ok time to do the stacking. Our first job is to do a nearest neighbour
309  stack in order to find the objects that will be rejected. If we
310  only want the NN algorithm, then we're done. Otherwise we can use the
311  rejection information to restack using any other algorithm */
312 
313  stack_nn();
314 
315  /* Switch for final stacking method */
316 
317  switch (method) {
318  case 0:
319  break;
320  case 1:
321  stack_lin();
322  break;
323  default:
324  break;
325  }
326 
327  /* Add the provenance cards to the header */
328 
329  vircam_prov(vircam_fits_get_ehu(outfits),inf,nimages);
330 
331  /* Assign some output info, tidy and get out of here */
332 
333  *out = outfits;
334  *outc = outconf;
335  tidy();
336  return(VIR_OK);
337 }
338 
339 static void seeing_wt(void) {
340  float *seeing;
341  int ierr,i;
342  cpl_propertylist *plist;
343  dstrct *dd;
344  const char *fctid = "seeing_wt";
345 
346  /* Get an array for the seeing estimates */
347 
348  seeing = cpl_malloc(nim*sizeof(*seeing));
349 
350  /* Read the seeing keyword. Get out of here if one is missing and
351  just use the background weights that we already have */
352 
353  ierr = 0;
354  for (i = 0; i < nim; i++) {
355  dd = fileptrs+i;
356  plist = vircam_fits_get_ehu(dd->fname);
357  if (cpl_propertylist_has(plist,"ESO DRS SEEING")) {
358  seeing[i] = cpl_propertylist_get_float(plist,"ESO DRS SEEING");
359  if (seeing[i] <= 0.0) {
360  ierr = 1;
361  break;
362  }
363  } else {
364  ierr = 1;
365  break;
366  }
367  }
368 
369  /* If there were no errors, then modify the background weights by the
370  seeing ratio */
371 
372  if (! ierr) {
373  sumweight = fileptrs->weight;
374  for (i = 1; i < nim; i++) {
375  dd = fileptrs + i;
376  dd->weight *= seeing[0]/seeing[i];
377  sumweight += dd->weight;
378  }
379  } else {
380  cpl_msg_warning(fctid,"No seeing weighting will be done");
381  }
382  freespace(seeing);
383 }
384 
385 /*---------------------------------------------------------------------------*/
401 /*---------------------------------------------------------------------------*/
402 
403 static int stack_lin(void) {
404  cpl_image *outim,*outimc;
405  double oxf,oyf;
406  int npts,i,*cdata,ix,iy,ind1,ind,ox,oy,m,jj,ii,nn,nx,ny,n;
407  float *workbufw,*data,wt,dx,dy,w[4],value,cval;
408  float med,noise;
409  dstrct *dd;
410 
411  /* Output file info */
412 
413  outim = vircam_fits_get_image(outfits);
414  outimc = vircam_fits_get_image(outconf);
415  nxo = (int)cpl_image_get_size_x(outim);
416  nyo = (int)cpl_image_get_size_y(outim);
417  npts = nxo*nyo;
418 
419  /* Get workspace for collecting info over the whole output map */
420 
421  workbufw = cpl_calloc(nxo*nyo,sizeof(*workbufw));
422 
423  /* Get output buffers */
424 
425  outdata = cpl_image_get_data_float(outim);
426  outdatac = cpl_image_get_data_int(outimc);
427 
428  /* Initialise the output map. This is because the output data array
429  probably still has the results of the nearest neighbour stack */
430 
431  for (i = 0; i < npts; i++) {
432  outdata[i] = 0.0;
433  outdatac[i] = 0;
434  }
435 
436  /* Ok, loop for each of the input files and collect all the information
437  for each pixel */
438 
439  for (i = 0; i < nim; i++) {
440  dd = fileptrs + i;
441  data = dd->data;
442  cdata = dd->cdata;
443  wt = dd->weight;
444  nx = dd->nx;
445  ny = dd->ny;
446 
447  /* Loop for each pixel in the data array and work out where it belongs
448  in the output array */
449 
450  n = 0;
451  for (iy = 0; iy < ny; iy++) {
452  ind1 = iy*nx;
453  for (ix = 0; ix < nx; ix++) {
454  ind = ind1 + ix;
455  outloc(dd->vwcs,(double)ix,(double)iy,outwcs,dd->trans,
456  &oxf,&oyf);
457  if (oxf < 0.0 || oyf < 0.0)
458  continue;
459  ox = (int)oxf;
460  oy = (int)oyf;
461  if ((dd->bpm)[ind] != 0)
462  continue;
463  if (ox < 0 || ox >= nxo || oy < 0 || oy >= nyo ||
464  cdata[ind] == 0)
465  continue;
466  dx = oxf - ox;
467  dy = oyf - oy;
468  w[0] = (1.0-dx)*(1.0-dy);
469  w[1] = dx*(1.0-dy);
470  w[2] = dy*(1.0-dx);
471  w[3] = dx*dy;
472  m = 0;
473  value = data[ind]/(dd->expscale) + dd->skydiff;
474  cval = (float)cdata[ind];
475  for (jj = 0; jj <= 1; jj++) {
476  for (ii = 0; ii <= 1; ii++) {
477  nn = (oy+jj)*nxo + ox + ii;
478  if (nn >= npts || ox+ii >= nxo || oy+jj >= nyo)
479  continue;
480  outdata[nn] += w[m]*value*wt*cval;
481  workbufw[nn] += w[m]*wt*cval;
482  m++;
483  }
484  }
485  }
486  }
487  }
488 
489  /* Now do the final averaging */
490 
491  for (i = 0; i < npts; i++) {
492  if (workbufw[i] != 0.0)
493  outdata[i] /= workbufw[i];
494  else
495  outdata[i] = fileptrs->sky;
496  workbufw[i] /= sumweight;
497  }
498  vircam_qmedsig(workbufw,NULL,npts,3.0,2,0.0,65535.0,&med,&noise);
499  for (i = 0; i < npts; i++) {
500  workbufw[i] *= (100.0/med);
501  outdatac[i] = max(0,min(1000,vircam_nint(workbufw[i])));
502  }
503  freespace(workbufw);
504  return(VIR_OK);
505 }
506 
507 /*---------------------------------------------------------------------------*/
523 /*---------------------------------------------------------------------------*/
524 
525 static int stack_nn(void) {
526  float *workbuf,*workbufc,*workbufw,lclip,hclip,*data,wt,med,noise,outvalc;
527  float value,*dbuf,*wbuf,*cbuf,outval,avlev,avvar,*confwork;
528  double oxf,oyf;
529  int *iloc,i,npts,ix,iy,ind1,ind,ox,oy,nn,bufloc,ncontrib,*lbuf;
530  int j,*cdata,nx,ny,n,nz;
531  short int clipped_low,clipped_high;
532  dstrct *dd;
533  unsigned char *id,*clipmap,*cc,*cclast,*ibuf;
534  short int *nbuf;
535  cpl_image *outim,*outimc;
536 
537  /* Output file info */
538 
539  outim = vircam_fits_get_image(outfits);
540  outimc = vircam_fits_get_image(outconf);
541  nxo = (int)cpl_image_get_size_x(outim);
542  nyo = (int)cpl_image_get_size_y(outim);
543 
544  /* Get workspace for collecting info over the whole output map */
545 
546  nz = nim + 2;
547  workbuf = cpl_calloc(nz*nxo*nyo,sizeof(*workbuf));
548  workbufc = cpl_calloc(nz*nxo*nyo,sizeof(*workbufc));
549  workbufw = cpl_calloc(nz*nxo*nyo,sizeof(*workbufw));
550  id = cpl_calloc(nz*nxo*nyo,sizeof(*id));
551  iloc = cpl_calloc(nz*nxo*nyo,sizeof(*iloc));
552  nbuf = cpl_calloc(nxo*nyo,sizeof(*nbuf));
553  confwork = cpl_calloc(nxo*nyo,sizeof(*confwork));
554 
555  /* Workspace for marking output rows with clipped pixels */
556 
557  clipmap = cpl_calloc(2*nxo,sizeof(*clipmap));
558 
559  /* Buffers for averaging an individual pixel */
560 
561  dbuf = cpl_calloc(2*nz,sizeof(*dbuf));
562  wbuf = cpl_calloc(2*nz,sizeof(*wbuf));
563  cbuf = cpl_calloc(2*nz,sizeof(*cbuf));
564  ibuf = cpl_calloc(2*nz,sizeof(*ibuf));
565  lbuf = cpl_calloc(2*nz,sizeof(*lbuf));
566 
567  /* Get output buffers */
568 
569  outdata = cpl_image_get_data_float(outim);
570  outdatac = cpl_image_get_data_int(outimc);
571 
572  /* Rejection thresholds */
573 
574  lclip = fileptrs->sky - lsig*(fileptrs->noise);
575  hclip = fileptrs->sky + hsig*(fileptrs->noise);
576 
577  /* Ok, loop for each of the input files and collect all the information
578  for each pixel */
579 
580  for (i = 0; i < nim; i++) {
581  dd = fileptrs + i;
582  nx = dd->nx;
583  ny = dd->ny;
584  npts = nx*ny;
585  data = dd->data;
586  cdata = dd->cdata;
587  wt = dd->weight;
588 
589  /* Loop for each pixel in the data array and work out where it belongs
590  in the output array */
591 
592  n = 0;
593  for (iy = 0; iy < ny; iy++) {
594  ind1 = iy*nx;
595  for (ix = 0; ix < nx; ix++) {
596  ind = ind1 + ix;
597  outloc(dd->vwcs,(double)ix,(double)iy,outwcs,dd->trans,
598  &oxf,&oyf);
599  ox = vircam_nint(oxf) - 1;
600  oy = vircam_nint(oyf) - 1;
601  if (ox < 0 || ox >= nxo || oy < 0 || oy >= nyo ||
602  cdata[ind] == 0)
603  continue;
604  nn = oy*nxo + ox;
605  bufloc = oy*nxo*nz + ox*nz + nbuf[nn];
606  value = data[ind]/(dd->expscale) + dd->skydiff;
607  workbuf[bufloc] = value;
608  workbufw[bufloc] = wt*(float)cdata[ind];
609  workbufc[bufloc] = (float)cdata[ind];
610  id[bufloc] = i;
611  nbuf[nn] += 1;
612  iloc[bufloc] = ind;
613  }
614  }
615  }
616 
617 
618  /* Now loop through the output arrays and form an initial estimate */
619 
620  for (iy = 0; iy < nyo; iy++) {
621  ind1 = iy*nxo;
622  cc = clipmap + (iy % 2)*nxo;
623  cclast = clipmap + ((iy-1) % 2)*nxo;
624  for (ix = 0; ix < nxo; ix++) {
625  ind = ind1 + ix;
626  bufloc = iy*nxo*nz + ix*nz;
627  ncontrib = nbuf[ind];
628  if (ncontrib == 0) {
629  outdata[ind] = fileptrs->sky;
630  outdatac[ind] = 0;
631  continue;
632  }
633  cc[ix] = 0;
634 
635  /* Put some stuff in buffers for the averaging routine */
636 
637  for (i = 0; i < ncontrib; i++) {
638  dbuf[i] = workbuf[bufloc+i];
639  wbuf[i] = workbufw[bufloc+i];
640  cbuf[i] = workbufc[bufloc+i];
641  ibuf[i] = id[bufloc+i];
642  lbuf[i] = iloc[bufloc+i];
643  }
644 
645  /* Do the averages with the nominal clipping */
646 
647  do_averages(ncontrib,dbuf,wbuf,cbuf,ibuf,lclip,hclip,0.0,
648  &clipped_low,&clipped_high,&outval,&outvalc);
649 
650  /* Store these away */
651 
652  outdata[ind] = outval;
653  confwork[ind] = outvalc;
654  cc[ix] = (clipped_high >= 0);
655  if (clipped_low >= 0)
656  ((fileptrs+ibuf[clipped_low])->bpm)[lbuf[clipped_low]] = 2;
657  if (clipped_high >= 0)
658  ((fileptrs+ibuf[clipped_high])->bpm)[lbuf[clipped_high]] = 1;
659  }
660 
661  /* Look at the non-edge pixels and see if the local variation can
662  show whether any clipping that has been done is justified. NB:
663  we are looking at the _previous_ row here */
664 
665  if (iy < 2)
666  continue;
667  for (ix = 1; ix < nxo-1; ix++) {
668  if (! cclast[ix])
669  continue;
670  ind = (iy-1)*nxo + ix;
671  bufloc = (iy-1)*nxo*nz + ix*nz;
672  ncontrib = nbuf[ind];
673  if (ncontrib == 0)
674  continue;
675  avlev = 0.0;
676  for (i = -1; i <= 1; i++)
677  for (j = -1; j <= 1; j++)
678  avlev += outdata[ind + i*nxo + j];
679  avlev /= 9.0;
680  avvar = 0.0;
681  for (i = -1; i <= 1; i++)
682  for (j = -1; j <= 1; j++)
683  avvar += (float)fabs(avlev - outdata[ind + i*nxo + j]);
684  avvar /= 9.0;
685  if (avlev <= hclip && avvar <= 2.0*(fileptrs->sky))
686  continue;
687 
688  /* Put some stuff in buffers for the averaging routine */
689 
690  for (i = 0; i < ncontrib; i++) {
691  dbuf[i] = workbuf[bufloc+i];
692  wbuf[i] = workbufw[bufloc+i];
693  cbuf[i] = workbufc[bufloc+i];
694  ibuf[i] = id[bufloc+i];
695  lbuf[i] = iloc[bufloc+i];
696  (fileptrs+ibuf[i])->bpm[lbuf[i]] = 0;
697  }
698 
699  /* Redo the averages with the nominal clipping */
700 
701  do_averages(ncontrib,dbuf,wbuf,cbuf,ibuf,lclip,hclip,avvar,
702  &clipped_low,&clipped_high,&outval,&outvalc);
703 
704  /* Store these away */
705 
706  outdata[ind] = outval;
707  confwork[ind] = outvalc;
708  if (clipped_low >= 0)
709  ((fileptrs+ibuf[clipped_low])->bpm)[lbuf[clipped_low]] = 2;
710  if (clipped_high >= 0)
711  ((fileptrs+ibuf[clipped_high])->bpm)[lbuf[clipped_high]] = 1;
712  }
713  }
714 
715  /* Now renormalise the confidence map */
716 
717  npts = nxo*nyo;
718  vircam_qmedsig(confwork,NULL,npts,3.0,2,25.0,65535.0,&med,&noise);
719  for (i = 0; i < npts; i++) {
720  confwork[i] *= (100.0/med);
721  outdatac[i] = max(0,min(1000,vircam_nint(confwork[i])));
722  }
723 
724  /* Free up some workspace */
725 
726  freespace(workbuf);
727  freespace(workbufc);
728  freespace(workbufw);
729  freespace(confwork);
730  freespace(id);
731  freespace(iloc);
732  freespace(nbuf);
733  freespace(clipmap);
734  freespace(dbuf);
735  freespace(wbuf);
736  freespace(cbuf);
737  freespace(ibuf);
738  freespace(lbuf);
739  return(VIR_OK);
740 }
741 
742 /*---------------------------------------------------------------------------*/
780 /*---------------------------------------------------------------------------*/
781 
782 static void do_averages(int ncontrib, float *data, float *wconf, float *conf,
783  unsigned char *id, float lclip, float hclip,
784  float extra, short int *lowcl, short int *highcl,
785  float *outval, float *outvalc) {
786  float minval,maxval,sum,wsum,csum,reflev,clipval,clipscale;
787  int imax,imin,i;
788 
789  /* Do the summation and keep track of the minimum and
790  maximum values */
791 
792  minval = 1.0e10;
793  maxval = -1.0e10;
794  sum = 0.0;
795  wsum = 0.0;
796  csum = 0.0;
797  imin = -1;
798  imax = -1;
799  for (i = 0; i < ncontrib; i++) {
800  sum += data[i]*wconf[i];
801  wsum += wconf[i];
802  csum += conf[i];
803  if (data[i] > maxval) {
804  maxval = data[i];
805  imax = i;
806  }
807  if (data[i] < minval) {
808  minval = data[i];
809  imin = i;
810  }
811  }
812 
813  /* First crack at output value */
814 
815  if (wsum > 0.0) {
816  *outval = sum/wsum;
817  } else {
818  *outval = fileptrs->sky;
819  }
820 
821  /* Look at high clipping for of cosmic rays */
822 
823  *highcl = -1;
824  if (maxval > hclip && wsum > 150.0 && csum > 150.0 &&
825  imin != -1 && imax != -1) {
826  reflev = (sum - data[imax]*wconf[imax])/(wsum - wconf[imax]);
827  clipscale = max(1.0,reflev)/max(1.0,(fileptrs+id[imax])->sky);
828  clipval = reflev + clipscale*hsig*(fileptrs+id[imax])->noise;
829  clipval += extra;
830 
831  /* Clip the maximum point */
832 
833  if (maxval > clipval) {
834  sum -= data[imax]*wconf[imax];
835  wsum -= wconf[imax];
836  *outval = reflev;
837  *highcl = imax;
838  }
839  }
840 
841  /* Clip low points (presumably serious cosmetic issues on detectors) */
842 
843  *lowcl = -1;
844  if (minval < lclip && wsum > 150.0 && csum > 150.0 &&
845  imin != -1 && imax != -1) {
846  reflev = (sum - data[imin]*wconf[imin])/(wsum - wconf[imin]);
847  clipscale = max(1.0,reflev)/max(1.0,(fileptrs+id[imin])->sky);
848  clipval = reflev + clipscale*hsig*(fileptrs+id[imin])->noise;
849  if (minval < clipval) {
850  sum -= data[imin]*wconf[imin];
851  wsum -= wconf[imin];
852  *outval = reflev;
853  *lowcl = imin;
854  }
855  }
856 
857  /* Output confidence (not normalised) */
858 
859  *outvalc = wsum/sumweight;
860 }
861 
862 /*---------------------------------------------------------------------------*/
876 /*---------------------------------------------------------------------------*/
877 
878 static void output_files(void) {
879  cpl_image *newim,*newconf;
880  double lowerleft[2],upperleft[2],lowerright[2],upperright[2],xout,yout;
881  double xmin,ymin,sh[2];
882  dstrct *dd;
883  int i,ixo,iyo;
884 
885  /* Initialise the first element of the results arrays since we're
886  doing all this relative to the first image */
887 
888  lowerleft[0] = 1.0;
889  lowerleft[1] = 1.0;
890  upperleft[0] = 1.0;
891  upperleft[1] = (double)fileptrs->ny;
892  lowerright[0] = (double)fileptrs->nx;
893  lowerright[1] = 1.0;
894  upperright[0] = (double)fileptrs->nx;
895  upperright[1] = (double)fileptrs->ny;
896 
897  /* Loop for all images relative to the first one */
898 
899  for (i = 1; i < nim; i++) {
900  dd = fileptrs + i;
901  outloc(dd->vwcs,(double)1.0,(double)1.0,fileptrs->vwcs,dd->trans,
902  &xout,&yout);
903  lowerleft[0] = min(lowerleft[0],xout);
904  lowerleft[1] = min(lowerleft[1],yout);
905  outloc(dd->vwcs,(double)1.0,(double)dd->ny,fileptrs->vwcs,dd->trans,
906  &xout,&yout);
907  upperleft[0] = min(upperleft[0],xout);
908  upperleft[1] = max(upperleft[1],yout);
909  outloc(dd->vwcs,(double)dd->nx,(double)1.0,fileptrs->vwcs,dd->trans,
910  &xout,&yout);
911  lowerright[0] = max(lowerright[0],xout);
912  lowerright[1] = min(lowerright[1],yout);
913  outloc(dd->vwcs,(double)dd->nx,(double)dd->ny,fileptrs->vwcs,dd->trans,
914  &xout,&yout);
915  upperright[0] = max(upperright[0],xout);
916  upperright[1] = max(upperright[1],yout);
917  }
918 
919  /* Ok, what are the limits? */
920 
921  ixo = vircam_nint(max(lowerright[0]-lowerleft[0],upperright[0]-upperleft[0])) + 1;
922  iyo = vircam_nint(max(upperright[1]-lowerright[1],upperleft[1]-lowerleft[1])) + 1;
923  xmin = min(lowerleft[0],upperleft[0]);
924  ymin = min(lowerleft[1],lowerright[1]);
925 
926  /* Create the image */
927 
928  newim = cpl_image_new((cpl_size)ixo,(cpl_size)iyo,CPL_TYPE_FLOAT);
929  outfits = vircam_fits_wrap(newim,fileptrs->fname,NULL,NULL);
930  cpl_propertylist_update_int(vircam_fits_get_ehu(outfits),"NAXIS1",ixo);
931  cpl_propertylist_update_int(vircam_fits_get_ehu(outfits),"NAXIS2",iyo);
932 
933  /* Update the reference point for the WCS */
934 
935  sh[0] = xmin;
936  sh[1] = ymin;
937  (void)vircam_crpixshift(vircam_fits_get_ehu(outfits),1.0,sh);
938  outwcs = getwcsinfo(vircam_fits_get_ehu(outfits));
939 
940  /* Now create the confidence map image */
941 
942  newconf = cpl_image_new((cpl_size)ixo,(cpl_size)iyo,CPL_TYPE_INT);
943  outconf = vircam_fits_wrap(newconf,outfits,NULL,NULL);
944 }
945 
946 /*---------------------------------------------------------------------------*/
969 /*---------------------------------------------------------------------------*/
970 
971 static void diffxy(mywcs *wref, mywcs *wprog, float *dx, float *dy) {
972  double xref,yref,xprog,yprog;
973 
974  /* Find the middle of the reference frame */
975 
976  xref = 0.5*(double)wref->nx;
977  yref = 0.5*(double)wref->ny;
978 
979  /* Work out the output position */
980 
981  outloc(wref,xref,yref,wprog,NULL,&xprog,&yprog);
982 
983  /* Now the offsets */
984 
985  *dx = (float)(xref - xprog);
986  *dy = (float)(yref - yprog);
987 }
988 
989 /*---------------------------------------------------------------------------*/
1018 /*---------------------------------------------------------------------------*/
1019 
1020 static void outloc(mywcs *win, double xin, double yin, mywcs *wout,
1021  cpl_array *trans, double *xout, double *yout) {
1022  double xt,yt,xi,eta,r,rfac,aa,tandec,denom,rp,*tdata;
1023  int i;
1024 
1025  /* Do the conversion. First to standard coordinates in the frame of
1026  the input image */
1027 
1028  xt = xin - win->crpix[0];
1029  yt = yin - win->crpix[1];
1030  xi = win->cd[0]*xt + win->cd[1]*yt;
1031  eta = win->cd[2]*xt + win->cd[3]*yt;
1032  if (fabs(xt) < 1.0e-6 && fabs(yt) < 1.0e-6) {
1033  rfac = 1.0;
1034  } else {
1035  rp = sqrt(xi*xi + eta*eta);
1036  r = rp;
1037  for (i = 0; i < 3; i++) {
1038  rfac = win->pv21 + win->pv23*pow(r,2.0) + win->pv25*pow(r,4.0);
1039  r = rp/rfac;
1040  }
1041  rfac = tan(r)/rp;
1042  }
1043  xi *= rfac;
1044  eta *= rfac;
1045  aa = atan(xi*win->secd/(1.0-eta*win->tand));
1046  if (xi != 0.0)
1047  tandec = (eta+win->tand)*sin(aa)/(xi*win->secd);
1048  else
1049  tandec = (eta+win->tand)/(1.0 - eta*win->tand);
1050 
1051  /* Now from standard coordinates in the frame of the output image */
1052 
1053  aa += (win->crval[0] - wout->crval[0]);
1054  denom = wout->tand*tandec + cos(aa);
1055  xi = wout->secd*sin(aa)/denom;
1056  eta = (tandec - wout->tand*cos(aa))/denom;
1057  rp = sqrt(xi*xi + eta*eta);
1058  r = atan(rp);
1059  rfac = wout->pv21 + wout->pv23*pow(rp,2.0) + wout->pv25*pow(rp,4.0);
1060  rfac *= r/rp;
1061  xi *= rfac;
1062  eta *= rfac;
1063  denom = wout->cd[0]*wout->cd[3] - wout->cd[1]*wout->cd[2];
1064  xt = (xi*wout->cd[3] - eta*wout->cd[1])/denom + wout->crpix[0];
1065  yt = (eta*wout->cd[0] - xi*wout->cd[2])/denom + wout->crpix[1];
1066 
1067  /* If a fine tune transformation is defined, then do that now */
1068 
1069  if (trans != NULL) {
1070  tdata = cpl_array_get_data_double(trans);
1071  *xout = xt*tdata[0] + yt*tdata[1] + tdata[2];
1072  *yout = xt*tdata[3] + yt*tdata[4] + tdata[5];
1073  } else {
1074  *xout = xt;
1075  *yout = yt;
1076  }
1077 }
1078 
1079 /*---------------------------------------------------------------------------*/
1095 /*---------------------------------------------------------------------------*/
1096 
1097 static void backgrnd_ov(void) {
1098  int i,dx,dy,*cdata,n,jj,ind1,nx,ii,ind;
1099  long npts1,npts,start[2],end[2];
1100  dstrct *dd1,*dd2;
1101  float *data,sky1,skynoise1,sky2,skynoise2,frac;
1102 
1103  /* Initialise a few things */
1104 
1105  dd1 = fileptrs;
1106  npts1 = dd1->nx*dd1->ny;
1107 
1108  /* Start by working out the full background of the first frame */
1109 
1110  skyest(dd1->data,dd1->cdata,npts1,3.0,&sky1,&skynoise1);
1111  dd1->sky = sky1;
1112  dd1->noise = skynoise1;
1113  dd1->sky /= dd1->expscale;
1114  dd1->skydiff = 0.0;
1115  dd1->weight = 1.0;
1116  sumweight = 1.0;
1117 
1118  /* Now loop for all the others */
1119 
1120  for (i = 1; i < nim; i++) {
1121  dd2 = fileptrs + i;
1122 
1123  /* Offset differences */
1124 
1125  dx = vircam_nint(dd1->xoff - dd2->xoff);
1126  dy = vircam_nint(dd1->yoff - dd2->yoff);
1127 
1128  /* Ok, here's the tricky bit. Work out the range of x and y for
1129  each of the images. The shift moves an image to the left and
1130  and down. Hence positive values of dx or dy moves the second
1131  frame further to the left of down compared to the first frame.
1132  Here are the coordinate ranges for the first image */
1133 
1134  start[0] = max(0,-dx);
1135  start[1] = max(0,-dy);
1136  end[0] = min(dd1->nx-1,(dd1->nx)-dx-1);
1137  end[1] = min(dd1->ny-1,(dd1->ny)-dy-1);
1138  nx = dd1->nx;
1139 
1140  /* How much does this cover the first image? */
1141 
1142  npts = (end[0]-start[0]+1)*(end[1]-start[1]+1);
1143  frac = (float)npts/(float)npts1;
1144 
1145  /* If the coverage is more than 50% then calculate the background
1146  in the coverage region for the first image */
1147 
1148  if (frac >= 0.5) {
1149  data = cpl_malloc(npts*sizeof(*data));
1150  cdata = cpl_malloc(npts*sizeof(*cdata));
1151  n = 0;
1152  for (jj = start[1]; jj <= end[1]; jj++) {
1153  ind1 = jj*nx;
1154  for (ii = start[0]; ii <= end[0]; ii++) {
1155  ind = ind1 + ii;
1156  data[n] = (dd1->data)[ind];
1157  cdata[n++] = (dd1->cdata)[ind];
1158  }
1159  }
1160  skyest(data,cdata,npts,3.0,&sky1,&skynoise1);
1161  freespace(data);
1162  freespace(cdata);
1163  } else {
1164  sky1 = dd1->sky;
1165  }
1166 
1167  /* And here are the coordinate ranges for the second image if
1168  the coverage is more than 50%. Get the subset you need. */
1169 
1170  if (frac > 0.5) {
1171  start[0] = max(0,dx);
1172  start[1] = max(0,dy);
1173  end[0] = min(dd2->nx-1,(dd2->nx)+dx-1);
1174  end[1] = min(dd2->ny-1,(dd2->ny)+dy-1);
1175  nx = dd2->nx;
1176  npts = (end[0]-start[0]+1)*(end[1]-start[1]+1);
1177  data = cpl_malloc(npts*sizeof(*data));
1178  cdata = cpl_malloc(npts*sizeof(*cdata));
1179  n = 0;
1180  for (jj = start[1]; jj <= end[1]; jj++) {
1181  ind1 = jj*nx;
1182  for (ii = start[0]; ii <= end[0]; ii++) {
1183  ind = ind1 + ii;
1184  data[n] = dd2->data[ind];
1185  cdata[n++] = dd2->cdata[ind];
1186  }
1187  }
1188  skyest(data,cdata,npts,3.0,&sky2,&skynoise2);
1189  freespace(data);
1190  freespace(cdata);
1191  } else {
1192  npts = (dd2->nx)*(dd2->ny);
1193  skyest(dd2->data,dd2->cdata,npts,3.0,&sky2,&skynoise2);
1194  }
1195 
1196  /* Store info away */
1197 
1198  dd2->sky = sky2/dd2->expscale;
1199  dd2->skydiff = sky1 - sky2;
1200  dd2->noise = skynoise2/dd2->expscale;
1201  if (dd2->noise != 0.0 && dd2->sky != 0.0)
1202  dd2->weight = (float)(pow(dd1->noise,2.0)/pow(dd2->noise,2.0));
1203  else
1204  dd2->weight = 0.0;
1205  sumweight += dd2->weight;
1206  }
1207 }
1208 
1209 /*---------------------------------------------------------------------------*/
1235 /*---------------------------------------------------------------------------*/
1236 
1237 static void skyest(float *data, int *cdata, long npts, float thresh,
1238  float *skymed, float *skynoise) {
1239  unsigned char *bpm;
1240  int i;
1241 
1242  /* Get a dummy bad pixel mask */
1243 
1244  bpm = cpl_calloc(npts,sizeof(*bpm));
1245  for (i = 0; i < npts; i++)
1246  bpm[i] = (cdata[i] == 0);
1247 
1248  /* Get the stats */
1249 
1250  vircam_qmedsig(data,bpm,npts,thresh,2,-1000.0,65535.0,skymed,skynoise);
1251 
1252  /* Clean up */
1253 
1254  freespace(bpm);
1255 }
1256 
1257 /*---------------------------------------------------------------------------*/
1279 /*---------------------------------------------------------------------------*/
1280 
1281 static cpl_table *mknewtab(vir_tfits *catref, dstrct *dd) {
1282  cpl_table *newtab;
1283  double xx,yy,xx2,yy2;
1284  int nr,i;
1285  float *x,*y;
1286 
1287  /* Make a copy of the input table */
1288 
1289  newtab = cpl_table_duplicate(vircam_tfits_get_table(catref));
1290  nr = (int)cpl_table_get_nrow(newtab);
1291 
1292  /* Get the columns of interest */
1293 
1294  x = cpl_table_get_data_float(newtab,"X_coordinate");
1295  y = cpl_table_get_data_float(newtab,"Y_coordinate");
1296 
1297  /* Loop for each object */
1298 
1299  for (i = 0; i < nr; i++) {
1300  xx = (double)x[i];
1301  yy = (double)y[i];
1302  outloc(fileptrs->vwcs,xx,yy,dd->vwcs,NULL,&xx2,&yy2);
1303  x[i] = (float)xx2;
1304  y[i] = (float)yy2;
1305  }
1306  return(newtab);
1307 }
1308 
1309 /*---------------------------------------------------------------------------*/
1324 /*---------------------------------------------------------------------------*/
1325 
1326 static cpl_array *transinit(void) {
1327  double *td;
1328  cpl_array *t;
1329 
1330  t = cpl_array_new(6,CPL_TYPE_DOUBLE);
1331  td = cpl_array_get_data_double(t);
1332  td[0] = 1.0;
1333  td[1] = 0.0;
1334  td[2] = 0.0;
1335  td[3] = 0.0;
1336  td[4] = 1.0;
1337  td[5] = 0.0;
1338  return(t);
1339 }
1340 
1341 /*---------------------------------------------------------------------------*/
1354 /*---------------------------------------------------------------------------*/
1355 
1356 static void fileptrs_close(void) {
1357  int i;
1358 
1359  for (i = 0; i < nim; i++) {
1360  freespace((fileptrs+i)->vwcs);
1361  freearray((fileptrs+i)->trans);
1362  }
1363  freespace(fileptrs);
1364 }
1365 
1366 /*---------------------------------------------------------------------------*/
1382 /*---------------------------------------------------------------------------*/
1383 
1384 static mywcs *getwcsinfo(cpl_propertylist *plist) {
1385  mywcs *w;
1386  int i;
1387 
1388  /* Copy over the relevant info */
1389 
1390  w = cpl_malloc(sizeof(mywcs));
1391  (void)vircam_pfits_get_crval1(plist,&(w->crval[0]));
1392  (void)vircam_pfits_get_crval2(plist,&(w->crval[1]));
1393  (void)vircam_pfits_get_crpix1(plist,&(w->crpix[0]));
1394  (void)vircam_pfits_get_crpix2(plist,&(w->crpix[1]));
1395  (void)vircam_pfits_get_cd11(plist,&(w->cd[0]));
1396  (void)vircam_pfits_get_cd12(plist,&(w->cd[1]));
1397  (void)vircam_pfits_get_cd21(plist,&(w->cd[2]));
1398  (void)vircam_pfits_get_cd22(plist,&(w->cd[3]));
1399  (void)vircam_pfits_get_pv21(plist,&(w->pv21));
1400  (void)vircam_pfits_get_pv23(plist,&(w->pv23));
1401  (void)vircam_pfits_get_pv25(plist,&(w->pv25));
1402  (void)vircam_pfits_get_naxis1(plist,&(w->nx));
1403  (void)vircam_pfits_get_naxis2(plist,&(w->ny));
1404 
1405  /* Change angles to radians */
1406 
1407  for (i = 0; i < 2; i++)
1408  w->crval[i] /= DEGRAD;
1409  for (i = 0; i < 4; i++)
1410  w->cd[i] /= DEGRAD;
1411 
1412  /* Add a couple of convenience values */
1413 
1414  w->tand = tan(w->crval[1]);
1415  w->secd = 1.0/cos(w->crval[1]);
1416 
1417  /* Get out of here */
1418 
1419  return(w);
1420 }
1421 
1422 static void tidy(void) {
1423  fileptrs_close();
1424  freespace(outwcs);
1425 }
1426 
1427 /*
1428 
1429 $Log: not supported by cvs2svn $
1430 Revision 1.15 2012/01/15 17:40:09 jim
1431 Minor modifications to take into accout the changes in cpl API for v6
1432 
1433 Revision 1.14 2010/09/09 12:11:09 jim
1434 Fixed problems with docs that make doxygen barf
1435 
1436 Revision 1.13 2010/06/03 12:15:31 jim
1437 A few mods to get rid of compiler warnings
1438 
1439 Revision 1.12 2009/06/25 14:14:50 jim
1440 Fixed a little error handling problem
1441 
1442 Revision 1.11 2009/06/01 04:57:19 jim
1443 tidied up a few things
1444 
1445 Revision 1.10 2009/05/31 05:54:55 jim
1446 Nothing
1447 
1448 Revision 1.9 2009/05/31 05:42:55 jim
1449 reweighting linear interpolation
1450 
1451 Revision 1.8 2009/05/30 20:38:41 jim
1452 fixed another normalisation error
1453 
1454 Revision 1.7 2009/05/30 19:41:04 jim
1455 Fixed another small bug in the confidence map scaling
1456 
1457 Revision 1.6 2009/05/29 15:39:49 jim
1458 changed the way confidence map is normalised
1459 
1460 Revision 1.5 2009/02/20 10:52:49 jim
1461 Removed some superfluous declarations and tidied up some things to get
1462 rid of compiler warnings
1463 
1464 Revision 1.4 2009/01/28 10:55:10 jim
1465 Fixed so that if the calculation of trans array fails, then it issues a
1466 warning, but keeps going
1467 
1468 Revision 1.3 2009/01/19 14:35:43 jim
1469 Modified search radii in matching routines
1470 
1471 Revision 1.2 2008/11/25 06:19:20 jim
1472 Added routine prologues
1473 
1474 Revision 1.1 2008/11/21 10:19:16 jim
1475 New routine
1476 
1477 
1478 */