VIRCAM Pipeline  2.3.12
vircam_utils.c
1 /* $Id: vircam_utils.c,v 1.77 2013-10-15 16:54:41 jim Exp $
2  *
3  * This file is part of the VIRCAM Pipeline
4  * Copyright (C) 2015 European Southern Observatory
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:54:41 $
24  * $Revision: 1.77 $
25  * $Name: not supported by cvs2svn $
26  */
27 
28 /* Includes */
29 
30 #ifdef HAVE_CONFIG_H
31 #include <config.h>
32 #endif
33 
34 #define _XOPEN_SOURCE
35 #include <sys/time.h>
36 #include <time.h>
37 #include <libgen.h>
38 #include <string.h>
39 #include <unistd.h>
40 
41 #include <cpl.h>
42 #include <math.h>
43 
44 #include <casu_utils.h>
45 #include <casu_stats.h>
46 #include <casu_fits.h>
47 #include <catalogue/imcore.h>
48 
49 #include "vircam_utils.h"
50 #include "vircam_pfits.h"
51 
52 
53 #define SZKEY 32
54 #define SZVAL 64
55 #define SZBUF 1024
56 #define INITALLOC 1024
57 #define NOMPIXSIZE 0.34
58 
59 typedef struct {
60  char *filter_name;
61  float atm_extcoef;
62  float mag_offset;
63  char **coleq_columns;
64  float *coleq_coefs;
65  float gal_extcoef;
66  float default_zp;
67  float default_zp_err;
68  int ncolumns_coleq;
69 } photstrct;
70 
71 /* Define some columns for illumination correction tables */
72 
73 #define NI_COLS 5
74 static const char *illcor_cols[NI_COLS] = {"xmin","xmax","ymin","ymax",
75  "illcor"};
76 
77 /* Static subroutine prototypes */
78 
79 static float madfunc(int npts, float *xt, float *yt, float b);
80 static double pixsize (cpl_propertylist *plist);
81 static int vircam_phot_open(cpl_table *phottab, char *filt, photstrct *p);
82 static void vircam_phot_close(photstrct *p);
83 static int extract_columns(cpl_table *tab, photstrct *p);
84 static int extract_coleq(cpl_table *tab, photstrct *p);
85 
99 /*---------------------------------------------------------------------------*/
114 /*---------------------------------------------------------------------------*/
115 
116 extern const char *vircam_get_license(void) {
117  const char *vircam_license =
118  "This file is part of the VIRCAM Instrument Pipeline\n"
119  "Copyright (C) 2015 European Southern Observatory\n"
120  "\n"
121  "This program is free software; you can redistribute it and/or modify\n"
122  "it under the terms of the GNU General Public License as published by\n"
123  "the Free Software Foundation; either version 2 of the License, or\n"
124  "(at your option) any later version.\n"
125  "\n"
126  "This program is distributed in the hope that it will be useful,\n"
127  "but WITHOUT ANY WARRANTY; without even the implied warranty of\n"
128  "MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the\n"
129  "GNU General Public License for more details.\n"
130  "\n"
131  "You should have received a copy of the GNU General Public License\n"
132  "along with this program; if not, write to the Free Software\n"
133  "Foundation, Inc., 59 Temple Place, Suite 330, Boston, \n"
134  "MA 02111-1307 USA";
135  return(vircam_license);
136 }
137 
138 /*---------------------------------------------------------------------------*/
163 /*---------------------------------------------------------------------------*/
164 
165 extern void vircam_exten_range(int inexten, const cpl_frame *fr, int *out1,
166  int *out2) {
167  int nvircam = 16,n,nmax;
168  const char *fctid = "vircam_exten_range";
169 
170  /* Right, how many frames actually exist? */
171 
172  n = cpl_frame_get_nextensions(fr);
173 
174  /* If there are less than the cannonical number, then issue a warning
175  message here */
176 
177  if (n < nvircam) {
178  cpl_msg_warning(fctid,
179  "Only %" CPL_SIZE_FORMAT " extensions out of %" CPL_SIZE_FORMAT " are present",
180  (cpl_size)n,(cpl_size)nvircam);
181  nmax = n;
182  } else {
183  nmax = nvircam;
184  }
185 
186  /* If the extension number requested is zero, then do all the extensions
187  that are available */
188 
189  if (inexten == 0) {
190  *out1 = 1;
191  *out2 = nmax;
192 
193  /* Otherwise, check that the requested extension is actually present */
194 
195  } else {
196 
197  if (inexten > nmax) {
198  cpl_msg_error(fctid,
199  "Requested extension %" CPL_SIZE_FORMAT " is not present",
200  (cpl_size)inexten);
201  *out1 = -1;
202  *out2 = -1;
203  } else {
204  *out1 = inexten;
205  *out2 = inexten;
206  }
207  }
208  return;
209 }
210 
211 /*---------------------------------------------------------------------------*/
237 /*---------------------------------------------------------------------------*/
238 
239 extern void vircam_madfit(int npts, float *xdata, float *ydata,
240  float *intercept, float *slope) {
241  int j;
242  float sx,sy,sxx,sxy,det,aa,bb,temp,chisq,sigb,b1,f1,b2,f2,f;
243 
244  /* Do a linear least squares for a first estimate */
245 
246  sx = 0.0;
247  sy = 0.0;
248  sxx = 0.0;
249  sxy = 0.0;
250  for (j = 0; j < npts; j++) {
251  sx += xdata[j];
252  sy += ydata[j];
253  sxx += xdata[j]*xdata[j];
254  sxy += xdata[j]*ydata[j];
255  }
256  det = (float)npts*sxx - sx*sx;
257  if (det == 0.0) {
258  *slope = 0.0;
259  *intercept = 0.0;
260  return;
261  }
262  aa = (sxx*sy - sx*sxy)/det;
263  bb = ((float)npts*sxy - sx*sy)/det;
264  chisq = 0.0;
265  for (j = 0; j < npts; j++) {
266  temp = ydata[j] - (aa + bb*xdata[j]);
267  chisq += temp*temp;
268  }
269  sigb = sqrt(chisq/det);
270  if (sigb == 0.0) {
271  *slope = bb;
272  *intercept = aa;
273  return;
274  }
275 
276  /* Now bracket the solution */
277 
278  b1 = bb;
279  f1 = madfunc(npts,xdata,ydata,b1);
280  b2 = bb + ((f1 > 0.0) ? fabs(3.0*sigb) : -fabs(3.0*sigb));
281  f2 = madfunc(npts,xdata,ydata,b2);
282  while (f1*f2 > 0.0) {
283  bb = 2.0*b2 - b1;
284  b1 = b2;
285  f1 = f2;
286  b2 = bb;
287  f2 = madfunc(npts,xdata,ydata,b2);
288  }
289 
290  /* Iterate to find the minimum value */
291 
292  sigb = 0.01*sigb;
293  while (fabs(b2 - b1) > sigb) {
294  bb = 0.5*(b1 + b2);
295  if (bb == b1 || bb == b2)
296  break;
297  f = madfunc(npts,xdata,ydata,bb);
298  if (f*f1 >= 0.0) {
299  f1 = f;
300  b1 = bb;
301  } else {
302  f2 = f;
303  b2 = bb;
304  }
305  }
306  *intercept = aa;
307  *slope = bb;
308 }
309 
310 /*---------------------------------------------------------------------------*/
333 /*---------------------------------------------------------------------------*/
334 
335 
336 static float madfunc(int npts, float *xt, float *yt, float b) {
337  float *arr,aa,d,sum;
338  int j;
339 
340  arr = cpl_malloc(npts*sizeof(*arr));
341  for (j = 0; j < npts; j++)
342  arr[j] = yt[j] - b*xt[j];
343  aa = casu_med(arr,NULL,(long)npts);
344  sum = 0.0;
345  for (j = 0; j < npts; j++) {
346  d = yt[j] - (b*xt[j] + aa);
347  sum += d > 0.0 ? xt[j] : -xt[j];
348  }
349  cpl_free(arr);
350  return(sum);
351 }
352 
353 /*---------------------------------------------------------------------------*/
380 /*---------------------------------------------------------------------------*/
381 
382 extern void vircam_linfit(int npts, double *xdata, double *ydata,
383  double *intercept, double *slope, double *sig) {
384  int j;
385  double sx,sy,sxx,sxy,det,aa,bb,temp,sum,sumsq;
386 
387  /* Do a linear least squares */
388 
389  sx = 0.0;
390  sy = 0.0;
391  sxx = 0.0;
392  sxy = 0.0;
393  for (j = 0; j < npts; j++) {
394  sx += xdata[j];
395  sy += ydata[j];
396  sxx += xdata[j]*xdata[j];
397  sxy += xdata[j]*ydata[j];
398  }
399  det = (double)npts*sxx - sx*sx;
400  if (det == 0.0) {
401  *slope = 0.0;
402  *intercept = 0.0;
403  *sig = 0.0;
404  return;
405  }
406 
407  /* Intercept and slope */
408 
409  aa = (sxx*sy - sx*sxy)/det;
410  bb = ((double)npts*sxy - sx*sy)/det;
411 
412  /* Work out RMS of fit */
413 
414  sum = 0.0;
415  sumsq = 0.0;
416  for (j = 0; j < npts; j++) {
417  temp = ydata[j] - (aa + bb*xdata[j]);
418  sum += temp;
419  sumsq += temp*temp;
420  }
421  sum /= (double)npts;
422 
423  /* Set return values */
424 
425  *sig = sqrt(sumsq/(double)npts - sum*sum);
426  *slope = bb;
427  *intercept = aa;
428 }
429 
430 /*---------------------------------------------------------------------------*/
454 /*---------------------------------------------------------------------------*/
455 
456 extern int vircam_solve_gauss(double **a, double *b, int m) {
457  double temp,big,pivot,rmax;
458  int i,iu,j,k,jl,ib,ir;
459  int l = 0;
460 
461  iu = m - 1;
462  for (i = 0; i < iu; i++) {
463  big = 0.0;
464 
465  /* find largest remaining term in ith column for pivot */
466 
467  for (k = i; k < m; k++) {
468  rmax = fabs(a[i][k]);
469  if (rmax > big) {
470  big = rmax;
471  l = k;
472  }
473  }
474 
475  /* check for non-zero term */
476 
477  if (big == 0.0) {
478  for (ib = 0; ib < m; ib++)
479  b[ib] = 0.0;
480  cpl_msg_error("vircam_solve_gauss","Zero Determinant\n");
481  return(CASU_FATAL);
482  }
483 
484  if (i != l) {
485 
486  /* switch rows */
487 
488  for (j = 0; j < m; j++) {
489  temp = a[j][i];
490  a[j][i] = a[j][l];
491  a[j][l] = temp;
492  }
493  temp = b[i];
494  b[i] = b[l];
495  b[l] = temp;
496  }
497 
498 
499  /* pivotal reduction */
500 
501  pivot = a[i][i];
502  jl = i+1;
503 
504  for (j = jl; j < m; j++) {
505  temp = a[i][j]/pivot;
506  b[j] -= temp*b[i];
507  for (k = i; k < m; k++)
508  a[k][j] -= temp*a[k][i];
509  }
510  }
511 
512  /* back substitution for solution */
513 
514  for (i = 0; i < m; i++) {
515  ir = m - 1 - i;
516  if (a[ir][ir] != 0.0) {
517  temp = b[ir];
518  if (ir != m - 1) {
519  for (j = 1; j <= i; j++) {
520  k = m - j;
521  temp -= a[k][ir]*b[k];
522  }
523  }
524  b[ir] = temp/a[ir][ir];
525  } else
526  b[ir] = 0.0;
527  }
528  return(CASU_OK);
529 }
530 
531 /*---------------------------------------------------------------------------*/
569 /*---------------------------------------------------------------------------*/
570 
571 extern int vircam_polyfit(const cpl_array *xarray, const cpl_array *yarray,
572  int ncoefs, int ilim, int niter, float lclip,
573  float hclip, cpl_array **polycf, double *sigfit) {
574  const char *fctid = "vircam_polyfit";
575  int npts,iter,i,j,nnew,k,retval,n;
576  double *xdata,*ydata,*pdata,*res,**a,*b,temp,sum,sumsq,val;
577  double lcut,hcut;
578  unsigned char *pm;
579 
580  /* Initialise a few things */
581 
582  *polycf = NULL;
583  *sigfit = -1.0;
584 
585  /* How many data points do we have? */
586 
587  npts = (int)cpl_array_get_size(xarray);
588 
589  /* Do we have enough points for the required order of the fit? */
590 
591  if (npts < ncoefs) {
592  cpl_msg_warning(fctid,
593  "Not data for fit, Npts = %" CPL_SIZE_FORMAT ", Ncoefs = %" CPL_SIZE_FORMAT,
594  (cpl_size)npts,(cpl_size)ncoefs);
595  return(CASU_FATAL);
596  }
597 
598  /* Create some arrays */
599 
600  a = cpl_malloc(ncoefs*sizeof(double *));
601  b = cpl_calloc(ncoefs,sizeof(double));
602  for (i = 0; i < ncoefs; i++)
603  a[i] = cpl_calloc(ncoefs,sizeof(double));
604  pm = cpl_calloc(npts,sizeof(unsigned char));
605  res = cpl_malloc(npts*sizeof(double));
606 
607  /* Get pointers to the input arrays */
608 
609  xdata = (double *)cpl_array_get_data_double_const(xarray);
610  ydata = (double *)cpl_array_get_data_double_const(yarray);
611 
612  /* Get an output array */
613 
614  *polycf = cpl_array_new((cpl_size)ncoefs,CPL_TYPE_DOUBLE);
615  pdata = cpl_array_get_data_double(*polycf);
616 
617  /* Iteration loop */
618 
619  for (iter = 0; iter <= niter; iter++) {
620 
621  /* Zero some accumulators */
622 
623  for (i = 0; i < ncoefs; i++) {
624  for (j = 0; j < ncoefs; j++)
625  a[i][j] = 0.0;
626  b[i] = 0.0;
627  }
628  nnew = 0;
629 
630  /* Cumulate sums */
631 
632  for (i = 0; i < npts; i++) {
633  if (pm[i] == 1)
634  continue;
635  for (k = 0; k < ncoefs; k++) {
636  temp = 1.0;
637  if (k + ilim != 0)
638  temp = pow(xdata[i],(double)(k+ilim));
639  b[k] += ydata[i]*temp;
640  for (j = 0; j <= k; j++) {
641  temp = 1.0;
642  if (k + j + 2*ilim != 0)
643  temp = pow(xdata[i],(double)(k+j+2*ilim));
644  a[j][k] += temp;
645  }
646  }
647  }
648  for (k = 1; k < ncoefs; k++)
649  for (j = 0; j < k; j++)
650  a[k][j] = a[j][k];
651 
652  /* Solve linear equations */
653 
654  retval = vircam_solve_gauss(a,b,ncoefs);
655  if (retval != CASU_OK) {
656  cpl_msg_warning(fctid,"Fit failed");
657  freearray(*polycf);
658  freespace2(a,ncoefs);
659  freespace(b);
660  freespace(pm);
661  freespace(res);
662  return(CASU_FATAL);
663  }
664 
665  /* Ok, assuming this is OK, then fill the polynomial coefficients */
666 
667  for (i = 0; i < ncoefs; i++)
668  pdata[i] = b[i];
669 
670  /* Calculate the fit quality */
671 
672  sum = 0.0;
673  sumsq = 0.0;
674  n = 0;
675  for (i = 0; i < npts; i++) {
676  if (pm[i] == 1)
677  continue;
678  val = 0.0;
679  for (j = 0; j < ncoefs; j++)
680  val += pdata[j]*pow(xdata[i],(double)j+ilim);
681  res[i] = val - ydata[i];
682  sum += res[i];
683  sumsq += pow(res[i],2.0);
684  n++;
685  }
686  sum /= (double)n;
687  *sigfit = sqrt(sumsq/(double)n - sum*sum);
688 
689  /* If this is not the last iteration, then do some clipping */
690 
691  lcut = sum - lclip*(*sigfit);
692  hcut = sum + hclip*(*sigfit);
693  if (iter < niter) {
694  for (i = 0; i < npts; i++) {
695  if (pm[i] == 1)
696  continue;
697  if (res[i] > hcut || res[i] < lcut) {
698  nnew++;
699  pm[i] = 1;
700  }
701  }
702  }
703 
704  /* If no new points have been clipped, then get out of here now... */
705 
706  if (nnew == 0)
707  break;
708  }
709 
710  /* Tidy up and get out of here */
711 
712  freespace2(a,ncoefs);
713  freespace(b);
714  freespace(pm);
715  freespace(res);
716  return(CASU_OK);
717 }
718 
719 /*---------------------------------------------------------------------------*/
760 /*---------------------------------------------------------------------------*/
761 
762 extern void vircam_difference_image(cpl_image *master, cpl_image *prog,
763  unsigned char *bpm, cpl_table *chantab,
764  int ncells, int oper, float *global_diff,
765  float *global_rms, cpl_image **diffim,
766  cpl_table **diffimstats) {
767  float *ddata,*work,mean,sig,med,mad;
768  long nx,ny,npts;
769  int nrows,i,nc1,nc2,nr,ixmin,ixmax,iymin,iymax,cnum,cx,cy,idx,idy;
770  int icx,icy,indy1,indy2,indx1,indx2,jp,jcx,jj,jcy,ii,ncx,ncy;
771  const char *fctid = "vircam_difference_image";
772 
773  /* Initialise the output */
774 
775  *diffim = NULL;
776  *diffimstats = NULL;
777  *global_diff = 0.0;
778  *global_rms = 0.0;
779 
780  /* Is there a programme frame or a master? */
781 
782  if (prog == NULL || master == NULL)
783  return;
784 
785  /* Start by subtracting the master image from the programme image */
786 
787  switch (oper) {
788  case 1:
789  *diffim = cpl_image_subtract_create(prog,master);
790  break;
791  case 2:
792  *diffim = cpl_image_divide_create(prog,master);
793  break;
794  default:
795  *diffim = NULL;
796  cpl_msg_error(fctid,"Invalid operation requested %" CPL_SIZE_FORMAT,
797  (cpl_size)oper);
798  break;
799  }
800  if (*diffim == NULL)
801  return;
802 
803  /* Work out a median difference */
804 
805  ddata = cpl_image_get_data_float(*diffim);
806  nx = (int)cpl_image_get_size_x(*diffim);
807  ny = (int)cpl_image_get_size_y(*diffim);
808  npts = nx*ny;
809  casu_medmad(ddata,bpm,npts,global_diff,global_rms);
810  *global_rms *= 1.48;
811 
812  /* Is there a channel table? */
813 
814  if (chantab == NULL)
815  return;
816 
817  /* Before going any further, check that the channel table has all of
818  the columns we need */
819 
820  if (! cpl_table_has_column(chantab,"ixmin") ||
821  ! cpl_table_has_column(chantab,"ixmax") ||
822  ! cpl_table_has_column(chantab,"iymin") ||
823  ! cpl_table_has_column(chantab,"iymax") ||
824  ! cpl_table_has_column(chantab,"channum")) {
825  cpl_msg_error(fctid,"Channel table is missing one of the required columns");
826 
827  return;
828  }
829 
830  /* Work out how to divide the channels */
831 
832  switch (ncells) {
833  case 1:
834  nc1 = 1;
835  nc2 = 1;
836  break;
837  case 2:
838  nc1 = 2;
839  nc2 = 1;
840  break;
841  case 4:
842  nc1 = 4;
843  nc2 = 1;
844  break;
845  case 8:
846  nc1 = 8;
847  nc2 = 1;
848  break;
849  case 16:
850  nc1 = 16;
851  nc2 = 1;
852  break;
853  case 32:
854  nc1 = 16;
855  nc2 = 2;
856  break;
857  case 64:
858  nc1 = 32;
859  nc2 = 2;
860  break;
861  default:
862  nc1 = 32;
863  nc2 = 2;
864  break;
865  }
866 
867  /* Create a difference image stats table */
868 
869  nrows = (int)cpl_table_count_selected(chantab);
870  *diffimstats = vircam_create_diffimg_stats(nrows*nc1*nc2);
871 
872  /* Loop for each data channel now */
873 
874  nr = 0;
875  for (i = 0; i < nrows; i++) {
876  ixmin = cpl_table_get_int(chantab,"ixmin",(cpl_size)i,NULL);
877  ixmax = cpl_table_get_int(chantab,"ixmax",(cpl_size)i,NULL);
878  iymin = cpl_table_get_int(chantab,"iymin",(cpl_size)i,NULL);
879  iymax = cpl_table_get_int(chantab,"iymax",(cpl_size)i,NULL);
880  cnum = cpl_table_get_int(chantab,"channum",(cpl_size)i,NULL);
881 
882  /* Which is the long axis? If the channels are rectangular then
883  divide the long axis by the greater number of cells. If the number
884  of cells is less than 8, then just divide the long axis. If the
885  channel is square, then divide the X axis more finely */
886 
887  cx = ixmax - ixmin + 1;
888  cy = iymax - iymin + 1;
889  if (cx > cy) {
890  ncx = max(nc1,nc2);
891  ncy = min(nc1,nc2);
892  } else if (cx < cy) {
893  ncy = max(nc1,nc2);
894  ncx = min(nc1,nc2);
895  } else {
896  ncx = max(nc1,nc2);
897  ncy = min(nc1,nc2);
898  }
899 
900  /* How big is the cell? */
901 
902  idy = cy/ncy;
903  idx = cx/ncx;
904  work = cpl_malloc(idx*idy*sizeof(*work));
905 
906  /* Now loop for aach cell */
907 
908  for (icy = 0; icy < ncy; icy++) {
909  indy1 = idy*icy;
910  indy2 = min(iymax,indy1+idy-1);
911  for (icx = 0; icx < ncx; icx++) {
912  indx1 = idx*icx;
913  indx2 = min(ixmax,indx1+idx-1);
914  jp = 0;
915  for (jcy = indy1; jcy < indy2; jcy++) {
916  jj = jcy*nx;
917  for (jcx = indx1; jcx < indx2; jcx++) {
918  ii = jj + jcx;
919  if (bpm != NULL && bpm[ii] == 0)
920  work[jp++] = ddata[ii];
921  }
922  }
923  (void)casu_meansig(work,NULL,(long)jp,&mean,&sig);
924  (void)casu_medmad(work,NULL,(long)jp,&med,&mad);
925  cpl_table_set_int(*diffimstats,"xmin",(cpl_size)nr,indx1+1);
926  cpl_table_set_int(*diffimstats,"xmax",(cpl_size)nr,indx2+1);
927  cpl_table_set_int(*diffimstats,"ymin",(cpl_size)nr,indy1+1);
928  cpl_table_set_int(*diffimstats,"ymax",(cpl_size)nr,indy2+1);
929  cpl_table_set_int(*diffimstats,"chan",(cpl_size)nr,cnum);
930  cpl_table_set_float(*diffimstats,"mean",(cpl_size)nr,mean);
931  cpl_table_set_float(*diffimstats,"median",(cpl_size)nr,med);
932  cpl_table_set_float(*diffimstats,"variance",(cpl_size)nr,
933  (sig*sig));
934  cpl_table_set_float(*diffimstats,"mad",(cpl_size)(nr++),mad);
935  }
936  }
937  cpl_free(work);
938  }
939 }
940 
941 /*---------------------------------------------------------------------------*/
958 /*---------------------------------------------------------------------------*/
959 
960 extern cpl_table *vircam_create_diffimg_stats(int nrows) {
961  cpl_table *diffimstats;
962 
963  diffimstats = cpl_table_new((cpl_size)nrows);
964  cpl_table_new_column(diffimstats,"xmin",CPL_TYPE_INT);
965  cpl_table_set_column_unit(diffimstats,"xmin","pixels");
966  cpl_table_new_column(diffimstats,"xmax",CPL_TYPE_INT);
967  cpl_table_set_column_unit(diffimstats,"xmax","pixels");
968  cpl_table_new_column(diffimstats,"ymin",CPL_TYPE_INT);
969  cpl_table_set_column_unit(diffimstats,"ymin","pixels");
970  cpl_table_new_column(diffimstats,"ymax",CPL_TYPE_INT);
971  cpl_table_set_column_unit(diffimstats,"ymax","pixels");
972  cpl_table_new_column(diffimstats,"chan",CPL_TYPE_INT);
973  cpl_table_set_column_unit(diffimstats,"chan","pixels");
974  cpl_table_new_column(diffimstats,"mean",CPL_TYPE_FLOAT);
975  cpl_table_set_column_unit(diffimstats,"mean","ADU");
976  cpl_table_new_column(diffimstats,"median",CPL_TYPE_FLOAT);
977  cpl_table_set_column_unit(diffimstats,"median","ADU");
978  cpl_table_new_column(diffimstats,"variance",CPL_TYPE_FLOAT);
979  cpl_table_set_column_unit(diffimstats,"variance","ADU**2");
980  cpl_table_new_column(diffimstats,"mad",CPL_TYPE_FLOAT);
981  cpl_table_set_column_unit(diffimstats,"mad","ADU");
982  return(diffimstats);
983 }
984 
985 /*---------------------------------------------------------------------------*/
1004 /*---------------------------------------------------------------------------*/
1005 
1006 extern int *vircam_dummy_confidence(long n) {
1007  int *cdata,i;
1008 
1009  cdata = cpl_malloc(n*sizeof(*cdata));
1010  for (i = 0; i < n; i++)
1011  cdata[i] = 100;
1012  return(cdata);
1013 }
1014 
1015 /*---------------------------------------------------------------------------*/
1031 /*---------------------------------------------------------------------------*/
1032 
1033 extern int vircam_is_dummy(cpl_propertylist *p) {
1034 
1035  /* Check for silly input */
1036 
1037  if (p == NULL)
1038  return(0);
1039 
1040  /* Check the propertylist and return the result */
1041 
1042  return(cpl_propertylist_has(p,"ESO DRS IMADUMMY"));
1043 }
1044 
1045 /*---------------------------------------------------------------------------*/
1062 /*---------------------------------------------------------------------------*/
1063 
1064 extern cpl_table *vircam_illcor_newtab(int nrows) {
1065  cpl_table *illcor;
1066  int i;
1067 
1068  /* Create the table structure and fill in the columns */
1069 
1070  illcor = cpl_table_new((cpl_size)nrows);
1071  for (i = 0; i < NI_COLS; i++)
1072  cpl_table_new_column(illcor,illcor_cols[i],CPL_TYPE_FLOAT);
1073  return(illcor);
1074 }
1075 
1076 /*---------------------------------------------------------------------------*/
1100 /*---------------------------------------------------------------------------*/
1101 
1102 extern void vircam_timestamp(char *out, int n) {
1103  struct timeval tv;
1104  struct tm *tm;
1105  float sec;
1106  struct tm result;
1107  /* Get the Greenwich Mean Time */
1108 
1109  (void)gettimeofday(&tv,NULL);
1110 
1111  tm = gmtime_r(&(tv.tv_sec), &result);
1112  sec = (float)tm->tm_sec + 1.0e-6*(float)tv.tv_usec;
1113 
1114  /* Now format it */
1115 
1116  (void)snprintf(out,n,"%04d-%02d-%02dT%02d:%02d:%07.4f",1900+tm->tm_year,
1117  tm->tm_mon+1,tm->tm_mday,tm->tm_hour,tm->tm_min,sec);
1118 }
1119 
1120 /*---------------------------------------------------------------------------*/
1146 /*---------------------------------------------------------------------------*/
1147 
1148 extern int vircam_check_crval(cpl_propertylist *phu, cpl_propertylist *ehu) {
1149  double crval1,crval2;
1150  cpl_property *p;
1151 
1152  /* Get the values of CRVAL1,2 from the extension header. If it doesn't
1153  exist, then signal an error */
1154 
1155  if ((vircam_pfits_get_crval1(ehu,&crval1) != CASU_OK) ||
1156  (vircam_pfits_get_crval2(ehu,&crval2) != CASU_OK))
1157  return(CASU_FATAL);
1158 
1159  /* If the values are both zero, then grab the values from the primary
1160  and silently update the extension header. NB: have to do this silly
1161  song and dance where we create new properties, rename the old ones,
1162  insert the new ones in the correct place and then erase the old
1163  ones because CPL won't do the update if it thinks the original
1164  CRVALs are floats rather than a double (which it inevitably will,
1165  since this is where CRVAL has been set identically to zero...) */
1166 
1167  if (fabs(crval1) < 1.0e-6 && fabs(crval2) < 1.0e-6) {
1168  if ((vircam_pfits_get_ra(phu,&crval1) != CASU_OK) ||
1169  (vircam_pfits_get_dec(phu,&crval2) != CASU_OK))
1170  return(CASU_FATAL);
1171  p = cpl_propertylist_get_property(ehu,"CRVAL1");
1172  cpl_property_set_name(p,"OLDCR1");
1173  p = cpl_propertylist_get_property(ehu,"CRVAL2");
1174  cpl_property_set_name(p,"OLDCR2");
1175  cpl_propertylist_insert_after_double(ehu,"OLDCR2","CRVAL1",crval1);
1176  cpl_propertylist_insert_after_double(ehu,"CRVAL1","CRVAL2",crval2);
1177  cpl_propertylist_erase(ehu,"OLDCR1");
1178  cpl_propertylist_erase(ehu,"OLDCR2");
1179  }
1180 
1181  /* Get out of here */
1182 
1183  return(CASU_OK);
1184 }
1185 
1186 /*---------------------------------------------------------------------------*/
1205 /*---------------------------------------------------------------------------*/
1206 
1207 extern void vircam_copywcs(cpl_propertylist *p1, cpl_propertylist *p2) {
1208  (void)cpl_propertylist_erase_regexp(p2,CPL_WCS_REGEXP,0);
1209  cpl_propertylist_copy_property_regexp(p2,p1,CPL_WCS_REGEXP,0);
1210 }
1211 
1212 /*---------------------------------------------------------------------------*/
1256 /*---------------------------------------------------------------------------*/
1257 
1258 extern void vircam_cull(cpl_frameset *in, float lthr, float hthr,
1259  int ndit, int jst, int jfn, cpl_array **mins,
1260  cpl_array **maxs, cpl_array **aves,
1261  cpl_frameset **out, int *status) {
1262  int i,n,nj,j,live,gotit;
1263  double val,dndit,fmin,fmax,sum;
1264  cpl_propertylist *p;
1265  cpl_image *im;
1266  cpl_frame *fr;
1267 
1268  /* Initialise a few things */
1269 
1270  nj = jfn - jst + 1;
1271  *status = 0;
1272  *out = NULL;
1273  gotit = 0;
1274  dndit = (double)ndit;
1275  *mins = cpl_array_new(nj,CPL_TYPE_DOUBLE);
1276  *maxs = cpl_array_new(nj,CPL_TYPE_DOUBLE);
1277  *aves = cpl_array_new(nj,CPL_TYPE_DOUBLE);
1278 
1279  /* Loop through each extension to be included */
1280 
1281  n = cpl_frameset_get_size(in);
1282  for (j = jst; j <= jfn; j++) {
1283 
1284  /* Look at the current extension in the first frame of the
1285  frameset. If it's flagged as dead, then assume the detector
1286  is dead for all of them. If this is the only extension we
1287  look at, then status = 1 is passed back to say all
1288  images are ignored because of a dead detector */
1289 
1290  fr = cpl_frameset_get_position(in,0);
1291  p = cpl_propertylist_load(cpl_frame_get_filename(fr),j);
1292  vircam_pfits_get_detlive(p,&live);
1293  cpl_propertylist_delete(p);
1294  if (! live) {
1295  cpl_array_set(*mins,j-jst,0.0);
1296  cpl_array_set(*maxs,j-jst,0.0);
1297  cpl_array_set(*aves,j-jst,0.0);
1298  if (nj == 1)
1299  *status = 1;
1300  continue;
1301  }
1302 
1303  /* OK, we have a live detector. If this is the first good detector
1304  then get an array we sort the output frameset here */
1305 
1306  if (! gotit) {
1307  *out = cpl_frameset_new();
1308  gotit = 1;
1309  }
1310 
1311  /* Loop through all of the images for the current detector and
1312  get an estimate of the median level. Keep tabs on the min
1313  and max. We normalise to ndit = 1 */
1314 
1315  fmin = 1000000000.0;
1316  fmax = -1000000000.0;
1317  sum = 0.0;
1318  for (i = 0; i < n; i++) {
1319  fr = cpl_frameset_get_position(in,i);
1320  im = cpl_image_load(cpl_frame_get_filename(fr),CPL_TYPE_FLOAT,0,j);
1321  val = cpl_image_get_median(im);
1322  val /= dndit;
1323  cpl_image_delete(im);
1324  fmin = min(fmin,val);
1325  fmax = max(fmax,val);
1326  sum += val;
1327  if (gotit == 1) {
1328  if (val >= lthr && val <= hthr)
1329  cpl_frameset_insert(*out,cpl_frame_duplicate(fr));
1330  }
1331  }
1332  cpl_array_set(*mins,j-jst,fmin);
1333  cpl_array_set(*maxs,j-jst,fmax);
1334  cpl_array_set(*aves,j-jst,sum/(double)n);
1335 
1336  /* Ok we've looked at the thresholds now. Flag so that we don't
1337  do it again. Check the number of frames that made it
1338  through the threshold and set the status to say if all of
1339  them failed */
1340 
1341  if (gotit == 1) {
1342  gotit = 2;
1343  if (cpl_frameset_get_size(*out) == 0)
1344  *status = 2;
1345  }
1346  }
1347 
1348  /* Right, get out of here */
1349 
1350  return;
1351 }
1352 
1353 /*---------------------------------------------------------------------------*/
1404 /*---------------------------------------------------------------------------*/
1405 
1406 extern int vircam_illum(casu_fits **images, cpl_table **mstds,
1407  cpl_propertylist **pl, int nimages, char *filt,
1408  cpl_table *phottab, int nbsize, cpl_table **illcor,
1409  float *illcor_rms,int *status) {
1410  const char *fctid = "vircam_illum";
1411  float **stdmagptr,fracx,fracy,**results,cdfudge,saturate,apcor3,exptime;
1412  float airmass,*catcore3,*xx,*yy,extinct,cf,fluxmag3,refmag,dm3,med,mad;
1413  float xmin,xmax,ymin,ymax,*resall,medall,lcut,hcut,sig,cut,*good;
1414  int nx,ny,ifracx,ifracy,nbsizx,nbsizy,nbx,nby,*nres,*nall,i,j,ncat,k;
1415  int ix,iy,ind,nresall,nallocall,ngood;
1416  casu_fits *im;
1417  cpl_propertylist *ehu_im,*ehu_cat;
1418  cpl_table *stds,*cl;
1419  photstrct p;
1420 
1421  /* Inherited status */
1422 
1423  *illcor = NULL;
1424  *illcor_rms = 0.0;
1425  if (*status != CASU_OK)
1426  return(*status);
1427 
1428  /* Check for nonsense errors */
1429 
1430  if (nimages <= 0) {
1431  cpl_msg_error(fctid,"No images included in photometric calibration");
1432  FATAL_ERROR
1433  }
1434 
1435  /* Set up the structure that will give us the colour equations for
1436  this filter later on */
1437 
1438  if (vircam_phot_open(phottab,filt,&p) != CASU_OK)
1439  FATAL_ERROR
1440 
1441  /* Get a workspace to hold the pointers to the magnitude columns */
1442 
1443  stdmagptr = cpl_malloc(p.ncolumns_coleq*sizeof(float *));
1444 
1445  /* Check to see if nbsize is close to exact divisor */
1446 
1447  nx = (int)cpl_image_get_size_x(casu_fits_get_image(images[0]));
1448  ny = (int)cpl_image_get_size_y(casu_fits_get_image(images[0]));
1449  fracx = ((float)nx)/((float)nbsize);
1450  fracy = ((float)ny)/((float)nbsize);
1451  ifracx = (int)(fracx + 0.1);
1452  ifracy = (int)(fracy + 0.1);
1453  nbsizx = nx/ifracx;
1454  nbsizy = ny/ifracy;
1455  nbsize = max(casu_nint(0.9*nbsize),min(nbsize,min(nbsizx,nbsizy)));
1456  nbsize = min(nx,min(ny,nbsize)); /* trap for small maps */
1457 
1458  /* Divide the map into partitions */
1459 
1460  nbx = nx/nbsize;
1461  nby = ny/nbsize;
1462 
1463  /* Get some workspace to hold all the zeropoints for all the images
1464  in the input list. This is an initial allocation and more can be
1465  made available later if needed */
1466 
1467  results = cpl_malloc(nbx*nby*sizeof(float *));
1468  good = cpl_malloc(nbx*nby*sizeof(float));
1469  nres = cpl_calloc(nbx*nby,sizeof(int));
1470  nall = cpl_malloc(nbx*nby*sizeof(int));
1471  for (i = 0; i < nbx*nby; i++) {
1472  results[i] = cpl_malloc(INITALLOC*sizeof(float));
1473  nall[i] = INITALLOC;
1474  }
1475  resall = cpl_malloc(INITALLOC*sizeof(float));
1476  nresall = 0;
1477  nallocall = INITALLOC;
1478 
1479  /* Set up the output illumination correction table */
1480 
1481  *illcor = vircam_illcor_newtab(nbx*nby);
1482 
1483  /* Loop for the input images and catalogues. Create some shorthand
1484  variables and work out what the zeropoint fudge factor will be
1485  based on the sampling of the current frame and the nominal sampling
1486  of a VIRCAM image */
1487 
1488  for (i = 0; i < nimages; i++) {
1489  im = images[i];
1490  ehu_im = casu_fits_get_ehu(im);
1491  stds = mstds[i];
1492  if (stds == NULL)
1493  continue;
1494  ehu_cat = pl[i];
1495  cdfudge = 2.5*log10((double)pixsize(ehu_im)/NOMPIXSIZE);
1496 
1497  /* Now for the input catalogues. Start by getting some useful info
1498  from the header */
1499 
1500  saturate = cpl_propertylist_get_float(ehu_cat,"ESO QC SATURATION");
1501  apcor3 = cpl_propertylist_get_float(ehu_cat,"APCOR3");
1502  (void)vircam_pfits_get_exptime(casu_fits_get_phu(im),&exptime);
1503  (void)vircam_pfits_get_airmass(casu_fits_get_phu(im),&airmass);
1504  if (cpl_error_get_code() != CPL_ERROR_NONE) {
1505  cpl_msg_error(fctid,"Unable to get header info for %s",
1507  cpl_error_reset();
1508  continue;
1509  }
1510 
1511  /* Get objects that are not saturated and aren't too elliptical */
1512 
1513  cpl_table_select_all(stds);
1514  (void)cpl_table_and_selected_float(stds,"Ellipticity",CPL_LESS_THAN,
1515  0.5);
1516  (void)cpl_table_and_selected_float(stds,"Peak_height",CPL_LESS_THAN,
1517  saturate);
1518  if (cpl_error_get_code() != CPL_ERROR_NONE) {
1519  cpl_msg_error(fctid,"Unable select data from matched stds tab %s",
1521  cpl_error_reset();
1522  continue;
1523  }
1524 
1525  /* Now extract all those that have passed the test. If none are left
1526  then issue a warning and move on to the next one */
1527 
1528  cl = cpl_table_extract_selected(stds);
1529  ncat = (int)cpl_table_get_nrow(cl);
1530  if (ncat == 0) {
1531  cpl_msg_error(fctid,"No good standards available for %s",
1533  cpl_table_delete(cl);
1534  cpl_error_reset();
1535  continue;
1536  }
1537 
1538  /* Dereference some of the columns */
1539 
1540  catcore3 = cpl_table_get_data_float(cl,"Aper_flux_3");
1541  xx = cpl_table_get_data_float(cl,"X_coordinate");
1542  yy = cpl_table_get_data_float(cl,"Y_coordinate");
1543  for (j = 0; j < p.ncolumns_coleq; j++)
1544  stdmagptr[j] = cpl_table_get_data_float(cl,(p.coleq_columns)[j]);
1545 
1546  /* Loop for all the standards */
1547 
1548  extinct = p.atm_extcoef*(airmass - 1.0);
1549  for (j = 0; j < ncat; j++) {
1550 
1551  /* Do core magnitude calculation */
1552 
1553  cf = catcore3[j]/exptime;
1554  if (cf < 1.0)
1555  cf = 1.0;
1556  fluxmag3 = 2.5*log10((double)cf) + apcor3;
1557 
1558  /* Work out a reference magnitude */
1559 
1560  refmag = p.mag_offset;
1561  for (k = 0; k < p.ncolumns_coleq; k++)
1562  refmag += ((p.coleq_coefs)[k]*stdmagptr[k][j]);
1563 
1564  /* Work out zero point */
1565 
1566  dm3 = refmag + fluxmag3 + extinct + cdfudge;
1567 
1568  /* What cell does this belong to? */
1569 
1570  ix = (int)(xx[j]/(float)nbsize);
1571  iy = (int)(yy[j]/(float)nbsize);
1572  ind = iy*nbx + ix;
1573  if (nres[ind] == nall[ind] - 1) {
1574  results[ind] = cpl_realloc(results[ind],
1575  (nall[ind]+INITALLOC)*sizeof(float));
1576  nall[ind] += INITALLOC;
1577  }
1578  results[ind][nres[ind]] = dm3;
1579  nres[ind] += 1;
1580 
1581  /* Add this into the global solution too */
1582 
1583  if (nresall == nallocall) {
1584  resall = cpl_realloc(resall,(nallocall+INITALLOC)*sizeof(float));
1585  nallocall += INITALLOC;
1586  }
1587  resall[nresall++] = dm3;
1588  }
1589 
1590  /* Get rid of the table subset */
1591 
1592  freetable(cl);
1593  }
1594 
1595  /* What is the global zero point ? */
1596 
1597  if (nresall != 0) {
1598  (void)casu_medmad(resall,NULL,(long)nresall,&medall,&mad);
1599  cut = max(3.0*1.48*mad,0.1);
1600  lcut = medall - cut;
1601  hcut = medall + cut;
1602  (void)casu_medmadcut(resall,NULL,(long)nresall,lcut,hcut,
1603  &medall,&sig);
1604  sig *= 1.48;
1605  }
1606 
1607  /* Ok, what is the mean zeropoint for each cell. Do all of this stuff
1608  even if nresall == 0 so that we at least have a table with all
1609  zeros at the end of this */
1610 
1611  ngood = 0;
1612  for (i = 0; i < nbx*nby; i++) {
1613  if (nres[i] > 0) {
1614  (void)casu_medmad(results[i],NULL,(long)nres[i],&med,&mad);
1615  cut = max(3.0*1.48*mad,0.3);
1616  lcut = med - cut;
1617  hcut = med + cut;
1618  (void)casu_medmadcut(results[i],NULL,(long)nres[i],lcut,hcut,
1619  &med,&sig);
1620  sig *= 1.48;
1621  med -= medall;
1622  good[ngood++] = med;
1623  } else {
1624  med = -99.0;
1625  }
1626 
1627  /* Where are we in the map? */
1628 
1629  iy = i/nbx;
1630  ix = i - iy*nbx;
1631  xmin = ix*nbsize;
1632  xmax = xmin + nbsize - 1;
1633  if (ix == nbx)
1634  xmax = nx;
1635  ymin = iy*nbsize;
1636  ymax = ymin + nbsize - 1;
1637  if (iy == nby)
1638  ymax = ny;
1639 
1640  /* Write out the results into illumination table */
1641 
1642  cpl_table_set_float(*illcor,"xmin",(cpl_size)i,xmin);
1643  cpl_table_set_float(*illcor,"xmax",(cpl_size)i,xmax);
1644  cpl_table_set_float(*illcor,"ymin",(cpl_size)i,ymin);
1645  cpl_table_set_float(*illcor,"ymax",(cpl_size)i,ymax);
1646  cpl_table_set_float(*illcor,"illcor",(cpl_size)i,med);
1647  }
1648 
1649  /* Get the RMS of the map */
1650 
1651  if (ngood > 0)
1652  casu_medsig(good,NULL,(long)ngood,&med,illcor_rms);
1653  else
1654  *illcor_rms = 0.0;
1655 
1656  /* Delete some workspace */
1657 
1658  for (i = 0; i < nbx*nby; i++)
1659  freespace(results[i]);
1660  freespace(results);
1661  freespace(nres);
1662  freespace(nall);
1663  freespace(stdmagptr);
1664  freespace(resall);
1665  freespace(good);
1666  vircam_phot_close(&p);
1667 
1668  /* Set out the appropriate return status */
1669 
1670  if (nresall != 0)
1671  GOOD_STATUS
1672  else
1673  WARN_RETURN
1674 
1675 }
1676 
1677 /*---------------------------------------------------------------------------*/
1702 /*---------------------------------------------------------------------------*/
1703 
1704 static int vircam_phot_open(cpl_table *phottab, char *filt, photstrct *p) {
1705  const char *fctid = "vircam_phot_open";
1706  int ns,nerr,null,nr;
1707  cpl_table *subset;
1708  char **filts;
1709  const char *req_cols[8] = {"filter_name","atm_extcoef","mag_offset",
1710  "coleq_columns","coleq_coefs","gal_extcoef",
1711  "default_zp","default_zp_err"};
1712 
1713  /* Initialise a few things */
1714 
1715  p->coleq_coefs = NULL;
1716  p->coleq_columns = NULL;
1717  p->ncolumns_coleq = 0;
1718  p->mag_offset = 0.0;
1719 
1720  /* Check the table and make sure it has all the required columns */
1721 
1722  nerr = 0;
1723  for (ns = 0; ns < 8; ns++) {
1724  if (! cpl_table_has_column(phottab,req_cols[ns])) {
1725  cpl_msg_error(fctid,"Photometry table missing column %s",
1726  req_cols[ns]);
1727  nerr++;
1728  }
1729  }
1730  if (nerr > 0)
1731  return(CASU_FATAL);
1732 
1733  /* Search the table and find the rows that matches the filter. NB: we
1734  can't use cpl_table_and_selected_string because it matches by regular
1735  expression rather than exact matching! So just loop through and choose
1736  the first one that matches exactly. There shouldn't be more than one! */
1737 
1738  filts = cpl_table_get_data_string(phottab,"filter_name");
1739  nr = cpl_table_get_nrow(phottab);
1740  nerr = 1;
1741  for (ns = 0; ns < nr; ns++) {
1742  if (strncmp(filts[ns],filt,16) == 0) {
1743  nerr = 0;
1744  break;
1745  }
1746  }
1747  if (nerr == 1) {
1748  cpl_msg_error(fctid,"Unable to match photometry table to filter %s",
1749  filt);
1750  return(CASU_FATAL);
1751  }
1752  cpl_table_and_selected_window(phottab,ns,1);
1753 
1754  /* Read the information from the selected row */
1755 
1756  subset = cpl_table_extract_selected(phottab);
1757  p->filter_name = (char *)cpl_table_get_string(subset,"filter_name",0);
1758  p->atm_extcoef = cpl_table_get_float(subset,"atm_extcoef",0,&null);
1759  p->mag_offset = cpl_table_get_float(subset,"mag_offset",0,&null);
1760  p->gal_extcoef = cpl_table_get_float(subset,"gal_extcoef",0,&null);
1761  p->default_zp = cpl_table_get_float(subset,"default_zp",0,&null);
1762  p->default_zp_err = cpl_table_get_float(subset,"default_zp_err",0,&null);
1763  if (extract_columns(subset,p) != CASU_OK) {
1764  freetable(subset);
1765  return(CASU_FATAL);
1766  }
1767  if (extract_coleq(subset,p) != CASU_OK) {
1768  freetable(subset);
1769  return(CASU_FATAL);
1770  }
1771  freetable(subset);
1772 
1773  /* Get out of here */
1774 
1775  return(CASU_OK);
1776 }
1777 
1778 /*---------------------------------------------------------------------------*/
1794 /*---------------------------------------------------------------------------*/
1795 
1796 static void vircam_phot_close(photstrct *p) {
1797  int j;
1798 
1799  for (j = 0; j < p->ncolumns_coleq; j++)
1800  freespace((p->coleq_columns)[j]);
1801  freespace(p->coleq_columns);
1802  freespace(p->coleq_coefs);
1803 }
1804 
1805 /*---------------------------------------------------------------------------*/
1824 /*---------------------------------------------------------------------------*/
1825 
1826 static int extract_columns(cpl_table *tab, photstrct *p) {
1827  int nv,i,j;
1828  char *v,*w,**z;
1829 
1830  /* Get the relevant value from the table */
1831 
1832  v = cpl_strdup(cpl_table_get_string(tab,"coleq_columns",0));
1833 
1834  /* Count the number of commas in the string to see how many
1835  columns there are */
1836 
1837  nv = 1;
1838  j = strlen(v);
1839  for (i = 0; i < j; i++)
1840  if (v[i] == ',')
1841  nv++;
1842  p->ncolumns_coleq = nv;
1843 
1844  /* Now parse the string into the column names */
1845 
1846  z = cpl_malloc(nv*sizeof(char *));
1847  char *saveptr = NULL;
1848  for (i = 0; i < nv; i++) {
1849  if (i == 0)
1850  w = strtok_r(v,",",&saveptr);
1851  else
1852  w = strtok_r(NULL,",",&saveptr);
1853  z[i] = cpl_strdup(w);
1854  }
1855  p->coleq_columns = z;
1856 
1857  freespace(v);
1858  return(CASU_OK);
1859 }
1860 
1861 /*---------------------------------------------------------------------------*/
1880 /*---------------------------------------------------------------------------*/
1881 
1882 static int extract_coleq(cpl_table *tab, photstrct *p) {
1883  int nv,i,j;
1884  char *v,*w;
1885  float *z;
1886 
1887  /* Get the relevant value from the table */
1888 
1889  v = cpl_strdup(cpl_table_get_string(tab,"coleq_coefs",0));
1890 
1891  /* Count the number of commas in the string to see how many
1892  columns there are */
1893 
1894  nv = 1;
1895  j = strlen(v);
1896  for (i = 0; i < j; i++)
1897  if (v[i] == ',')
1898  nv++;
1899 
1900  /* Now parse the string into the colour equation values */
1901 
1902  z = cpl_malloc(nv*sizeof(float));
1903  char *saveptr1;
1904  for (i = 0; i < nv; i++) {
1905  if (i == 0)
1906  w = strtok_r(v,",", &saveptr1);
1907  else
1908  w = strtok_r(NULL,",", &saveptr1);
1909  z[i] = (float)atof(w);
1910  }
1911  p->coleq_coefs = z;
1912  freespace(v);
1913  return(CASU_OK);
1914 }
1915 
1916 /*---------------------------------------------------------------------------*/
1933 /*---------------------------------------------------------------------------*/
1934 
1935 static double pixsize (cpl_propertylist *plist) {
1936  cpl_wcs *wcs;
1937  double *cd,pix;
1938 
1939  wcs = cpl_wcs_new_from_propertylist(plist);
1940  cd = cpl_matrix_get_data((cpl_matrix *)cpl_wcs_get_cd(wcs));
1941  pix = 3600.0*sqrt(fabs(cd[0]*cd[3] - cd[1]*cd[2]));
1942  cpl_wcs_delete(wcs);
1943  return(pix);
1944 }
1945 
1948 /*
1949 
1950 $Log: not supported by cvs2svn $
1951 Revision 1.76 2012/01/27 12:25:10 jim
1952 Fixed some casts
1953 
1954 Revision 1.75 2012/01/16 12:32:18 jim
1955 A few more changes to fit in with cpl6
1956 
1957 Revision 1.74 2012/01/15 17:40:09 jim
1958 Minor modifications to take into accout the changes in cpl API for v6
1959 
1960 Revision 1.73 2010/12/09 13:23:56 jim
1961 Minor doc fix
1962 
1963 Revision 1.72 2010/09/09 12:11:09 jim
1964 Fixed problems with docs that make doxygen barf
1965 
1966 Revision 1.71 2010/06/30 12:42:00 jim
1967 A few fixes to stop compiler compaints
1968 
1969 Revision 1.70 2010/06/07 12:42:40 jim
1970 Modifications to get rid of compiler gripes
1971 
1972 Revision 1.69 2010/03/22 06:07:24 jim
1973 Fixed bug in vircam_overexp
1974 
1975 Revision 1.68 2010/03/22 06:04:14 jim
1976 vircam_overexp now works out the ensemble stats on all frames rather than
1977 just those that are between the threshold levels.
1978 
1979 Revision 1.67 2010/01/31 19:30:45 jim
1980 Fixed bug in where global MAD wasn't being multiplied by 1.48 to turn it
1981 into an global RMS in the difference image routine
1982 
1983 Revision 1.66 2009/09/21 11:57:27 jim
1984 Modified vircam_overexp to pass back ensemble stats.
1985 
1986 Revision 1.65 2009/06/08 08:07:13 jim
1987 Fixed bugs in vircam_polyfit, where clipping loop was being done one too few
1988 times and where the error estimate was being done incorrectly
1989 
1990 Revision 1.64 2009/05/20 12:21:52 jim
1991 modified _propertylist_merge to avoid clash of data types
1992 
1993 Revision 1.63 2009/02/23 10:45:13 jim
1994 Fixed vircam_backmap so that if there were no points in a bin the it does
1995 something sensible to the value
1996 
1997 Revision 1.62 2009/02/20 11:35:11 jim
1998 Fixed bug in vircam_polyfit which stopped the iteration scheme.
1999 
2000 Revision 1.61 2008/12/08 06:39:45 jim
2001 Added new routine vircam_check_crval
2002 
2003 Revision 1.60 2008/11/25 18:55:36 jim
2004 rewrote vircam_sort to be a bit more reliable
2005 
2006 Revision 1.59 2008/11/02 14:35:07 jim
2007 Fixed teensy problem with two error messages
2008 
2009 Revision 1.58 2008/10/13 08:16:57 jim
2010 Fixed call to cpl_msg_warning
2011 
2012 Revision 1.57 2008/10/01 04:58:30 jim
2013 Added vircam_frameset_fexists
2014 
2015 Revision 1.56 2008/07/10 13:04:03 jim
2016 Fixed some comments
2017 
2018 Revision 1.55 2008/05/06 12:15:02 jim
2019 Changed vircam_catpars to check that we can read the index file and to send
2020 out error messages if there are problems
2021 
2022 Revision 1.54 2008/05/06 08:38:29 jim
2023 Modified call to cpl_propertylist_get_ from float to double in
2024 vircam_gaincor_calc.
2025 
2026 Revision 1.53 2007/11/14 14:47:53 jim
2027 vircam_linfit now works only with doubles
2028 
2029 Revision 1.52 2007/11/14 10:46:07 jim
2030 Fixed bugs in vircam_polyfit
2031 
2032 Revision 1.51 2007/10/25 17:33:31 jim
2033 Modified to add vircam_polyfit and to remove lint warnings
2034 
2035 Revision 1.50 2007/10/19 09:25:10 jim
2036 Fixed problems with missing includes
2037 
2038 Revision 1.49 2007/10/19 06:55:06 jim
2039 Modifications made to use new method for directing the recipes to the
2040 standard catalogues using the sof
2041 
2042 Revision 1.48 2007/10/15 12:50:28 jim
2043 Modified for compatibility with cpl_4.0
2044 
2045 Revision 1.47 2007/07/09 13:21:06 jim
2046 Changed argument list to vircam_exten_range to include a cpl_frame
2047 
2048 Revision 1.46 2007/05/08 10:41:01 jim
2049 Added vircam_gaincor_calc
2050 
2051 Revision 1.45 2007/05/03 11:15:33 jim
2052 Fixed little problem with table wcs
2053 
2054 Revision 1.44 2007/05/02 09:14:48 jim
2055 Added vircam_findcol and vircam_rename_property
2056 
2057 Revision 1.43 2007/03/01 12:42:42 jim
2058 Modified slightly after code checking
2059 
2060 Revision 1.42 2007/02/14 12:53:34 jim
2061 Removed vircam_paf_print_header
2062 
2063 Revision 1.41 2007/02/07 10:12:15 jim
2064 Removed vircam_ndit_correct as it is now no longer necessary
2065 
2066 Revision 1.40 2007/01/15 12:36:27 jim
2067 Fixed bug in nditcor where factor was being divided by instead of multiplied
2068 
2069 Revision 1.39 2006/12/06 12:58:41 jim
2070 Fixed a small bug in vircam_backmap
2071 
2072 Revision 1.38 2006/12/01 14:08:33 jim
2073 added vircam_backmap
2074 
2075 Revision 1.37 2006/11/28 20:52:17 jim
2076 Added vircam_illcor_newtab
2077 
2078 Revision 1.36 2006/11/27 12:00:17 jim
2079 cosmetic changes
2080 
2081 Revision 1.35 2006/10/05 09:23:00 jim
2082 Small modifications to a couple of cpl calls to bring them into line with
2083 cpl v3.0
2084 
2085 Revision 1.34 2006/08/27 20:28:15 jim
2086 added vircam_is_dummy
2087 
2088 Revision 1.33 2006/08/21 09:09:29 jim
2089 Added routines vircam_create_diffimg_stats and vircam_dummy_property
2090 
2091 Revision 1.32 2006/07/11 14:53:58 jim
2092 Fixed vircam_ndit_correct so that images with bad status are not corrected
2093 
2094 Revision 1.31 2006/07/04 09:19:06 jim
2095 replaced all sprintf statements with snprintf
2096 
2097 Revision 1.30 2006/06/30 21:32:09 jim
2098 Fixed bug in vircam_prov so that sprintf is replaced by snprintf
2099 
2100 Revision 1.29 2006/06/20 18:59:51 jim
2101 Fixed small bug in vircam_ndit_correct
2102 
2103 Revision 1.28 2006/06/09 22:24:47 jim
2104 Added vircam_ndit_correct
2105 
2106 Revision 1.27 2006/06/09 11:26:26 jim
2107 Small changes to keep lint happy
2108 
2109 Revision 1.26 2006/06/06 13:07:54 jim
2110 Change the stats routine used in vircam_difference_image
2111 
2112 Revision 1.25 2006/05/30 13:31:47 jim
2113 Added vircam_timestamp
2114 
2115 Revision 1.24 2006/05/24 13:35:49 jim
2116 Added vircam_dummy_image and vircam_dummy_catalogue
2117 
2118 Revision 1.23 2006/05/04 15:16:33 jim
2119 Fixed memory bug in vircam_overexp
2120 
2121 Revision 1.22 2006/04/20 11:29:41 jim
2122 Added vircam_overexp
2123 
2124 Revision 1.21 2006/03/22 11:42:32 jim
2125 Moved some routines into vircam_mask
2126 
2127 Revision 1.20 2006/03/15 10:43:42 jim
2128 Fixed a few things
2129 
2130 Revision 1.19 2006/03/08 14:32:22 jim
2131 Lots of little modifications
2132 
2133 Revision 1.18 2006/03/01 10:31:29 jim
2134 Now uses new vir_fits objects
2135 
2136 Revision 1.17 2006/02/18 11:48:55 jim
2137 *** empty log message ***
2138 
2139 Revision 1.16 2006/01/23 10:30:50 jim
2140 Mainly documentation mods
2141 
2142 Revision 1.15 2005/12/14 22:17:33 jim
2143 Updated docs
2144 
2145 Revision 1.14 2005/12/01 16:26:03 jim
2146 Added vircam_dummy_confidence
2147 
2148 Revision 1.13 2005/11/25 15:33:22 jim
2149 Some code fixes to keep splint happy
2150 
2151 Revision 1.12 2005/11/25 09:56:15 jim
2152 Tidied up some more documentation
2153 
2154 Revision 1.11 2005/11/25 09:36:23 jim
2155 Added vircam_linfit
2156 
2157 Revision 1.10 2005/11/07 13:15:16 jim
2158 Fixed lots of bugs and added some error checking
2159 
2160 Revision 1.9 2005/11/03 15:16:28 jim
2161 Lots of changes mainly to strengthen error reporting
2162 
2163 Revision 1.8 2005/11/03 13:28:50 jim
2164 All sorts of changes to tighten up error handling
2165 
2166 Revision 1.7 2005/10/14 13:21:04 jim
2167 Added some new routines
2168 
2169 Revision 1.6 2005/09/20 15:07:46 jim
2170 Fixed a few bugs and added a few things
2171 
2172 Revision 1.5 2005/09/07 12:47:22 jim
2173 renamed a malloc and free to the cpl equivallent
2174 
2175 Revision 1.4 2005/08/10 09:55:05 jim
2176 Modified vircam_madfit so that if determinant is zero, then it sends a
2177 zero slope and intercept back
2178 
2179 Revision 1.3 2005/08/10 08:42:00 jim
2180 Modified vircam_madfit so that if the initial least squares fit gets a
2181 perfect result, then it just returns with those values rather than to try
2182 and do the median fit
2183 
2184 Revision 1.2 2005/08/09 15:30:00 jim
2185 vircam_exten_range had number of chips wrong!
2186 
2187 Revision 1.1.1.1 2005/08/05 08:29:09 jim
2188 Initial import
2189 
2190 
2191 */
2192 
cpl_image * casu_fits_get_image(casu_fits *p)
Definition: casu_fits.c:436
char * casu_fits_get_fullname(casu_fits *p)
Definition: casu_fits.c:680
cpl_propertylist * casu_fits_get_phu(casu_fits *p)
Definition: casu_fits.c:531
cpl_propertylist * casu_fits_get_ehu(casu_fits *p)
Definition: casu_fits.c:576
void casu_medmad(float *data, unsigned char *bpm, long np, float *med, float *mad)
Definition: casu_stats.c:347
float casu_med(float *data, unsigned char *bpm, long npts)
Definition: casu_stats.c:89
void casu_medsig(float *data, unsigned char *bpm, long np, float *med, float *sig)
Definition: casu_stats.c:672
void casu_medmadcut(float *data, unsigned char *bpm, long np, float lcut, float hcut, float *med, float *mad)
Definition: casu_stats.c:406
int casu_meansig(float *data, unsigned char *bpm, long npts, float *mean, float *sig)
Definition: casu_stats.c:595
int vircam_illum(casu_fits **images, cpl_table **mstds, cpl_propertylist **pl, int nimages, char *filt, cpl_table *phottab, int nbsize, cpl_table **illcor, float *illcor_rms, int *status)
Work out the illumination correction.
int vircam_pfits_get_ra(const cpl_propertylist *plist, double *ra)
Get the value of RA.
Definition: vircam_pfits.c:745
int vircam_pfits_get_dec(const cpl_propertylist *plist, double *dec)
Get the value of DEC.
Definition: vircam_pfits.c:761
int vircam_pfits_get_crval1(const cpl_propertylist *plist, double *crval1)
Get the value of crval1.
Definition: vircam_pfits.c:68
int vircam_pfits_get_crval2(const cpl_propertylist *plist, double *crval2)
Get the value of crval2.
Definition: vircam_pfits.c:102
int vircam_pfits_get_detlive(const cpl_propertylist *plist, int *detlive)
Get the value of DET_LIVE.
Definition: vircam_pfits.c:624
int vircam_pfits_get_exptime(const cpl_propertylist *plist, float *exptime)
Get the value of exposure time.
Definition: vircam_pfits.c:245
int vircam_pfits_get_airmass(const cpl_propertylist *plist, float *airmass)
Get the value of the airmass.
Definition: vircam_pfits.c:405
void vircam_timestamp(char *out, int n)
cpl_table * vircam_illcor_newtab(int nrows)
int vircam_solve_gauss(double **a, double *b, int m)
Definition: vircam_utils.c:456
void vircam_linfit(int npts, double *xdata, double *ydata, double *intercept, double *slope, double *sig)
Definition: vircam_utils.c:382
int * vircam_dummy_confidence(long n)
void vircam_cull(cpl_frameset *in, float lthr, float hthr, int ndit, int jst, int jfn, cpl_array **mins, cpl_array **maxs, cpl_array **aves, cpl_frameset **out, int *status)
int vircam_check_crval(cpl_propertylist *phu, cpl_propertylist *ehu)
int vircam_polyfit(const cpl_array *xarray, const cpl_array *yarray, int ncoefs, int ilim, int niter, float lclip, float hclip, cpl_array **polycf, double *sigfit)
Definition: vircam_utils.c:571
cpl_table * vircam_create_diffimg_stats(int nrows)
Definition: vircam_utils.c:960
void vircam_copywcs(cpl_propertylist *p1, cpl_propertylist *p2)
void vircam_difference_image(cpl_image *master, cpl_image *prog, unsigned char *bpm, cpl_table *chantab, int ncells, int oper, float *global_diff, float *global_rms, cpl_image **diffim, cpl_table **diffimstats)
Definition: vircam_utils.c:762
void vircam_madfit(int npts, float *xdata, float *ydata, float *intercept, float *slope)
Definition: vircam_utils.c:239
int vircam_is_dummy(cpl_propertylist *p)
const char * vircam_get_license(void)
Definition: vircam_utils.c:116
void vircam_exten_range(int inexten, const cpl_frame *fr, int *out1, int *out2)
Definition: vircam_utils.c:165