VIRCAM Pipeline 2.3.15
casu_imstack.c
1/* $Id: casu_imstack.c,v 1.10 2015/11/25 10:27:19 jim Exp $
2 *
3 * This file is part of the CASU Pipeline utilities
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: 2015/11/25 10:27:19 $
24 * $Revision: 1.10 $
25 * $Name: $
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 <string.h>
37#include <unistd.h>
38#include "casu_mods.h"
39#include "catalogue/casu_utils.h"
40#include "catalogue/casu_fits.h"
41#include "casu_stats.h"
42#include "catalogue/casu_wcsutils.h"
43#include "catalogue/casu_tfits.h"
44
45#define DATAMIN -1000.0
46#define DATAMAX 65535.0
47
48typedef struct {
49 int nx;
50 int ny;
51 casu_fits *fname;
52 casu_fits *conf;
53 casu_fits *var;
54 unsigned char *bpm;
55 cpl_wcs *vwcs;
56 float xoff;
57 float yoff;
58 cpl_array *trans;
59 cpl_array *invtrans;
60 float sky;
61 float skydiff;
62 float noise;
63 float expscale;
64 float weight;
65 double ramin;
66 double ramax;
67 double decmin;
68 double decmax;
69 int confcopy;
70} dstrct;
71
72static void linecoverage(cpl_wcs *outwcs, cpl_matrix *xyin, double *lra1,
73 double *lra2, double *ldec1, double *ldec2);
74static int stack_nn(int nim, dstrct *fileptrs, float lsig, float hsig,
75 float sumweights, casu_fits *outfits, casu_fits *outconf,
76 casu_fits *outvar, cpl_wcs *outwcs);
77static int stack_lin(int nim, dstrct *fileptrs, float sumweights,
78 casu_fits *outfits, casu_fits *outconf, casu_fits *outvar,
79 cpl_wcs *outwcs);
80static void stackslow(int nimages, dstrct *fileptrs, float lsig, float hsig,
81 float sumweights, casu_fits *outfits, casu_fits *outconf,
82 casu_fits *outvar, cpl_wcs *outwcs, int dolin);
83static void seeing_wt(int nim, dstrct *fileptrs);
84static void output_files(int nim, dstrct *fileptrs, casu_fits **outfits,
85 casu_fits **outconf, casu_fits **outvar,
86 cpl_wcs **outwcs);
87static void diffxy(cpl_wcs *wref, cpl_wcs *wprog, float *dx, float *dy);
88static void outloc(cpl_wcs *win, cpl_matrix *in, cpl_wcs *wout,
89 cpl_array *trans, cpl_matrix **out);
90static void inloc(cpl_wcs *wout, cpl_matrix *xyout, cpl_wcs *win,
91 cpl_array *trans, cpl_matrix **xyin);
92static void backgrnd_ov(int nim, dstrct *fileptrs, int unmap);
93static void skyest(float *data, int *cdata, long npts, float thresh,
94 float *skymed, float *skynoise);
95static void do_averages(int ncontrib, float *data, float *vdata, float *wconf,
96 float *conf, unsigned char *id, float lclip,
97 float hclip, float hsig, float extra, dstrct *fileptrs,
98 float sumweights, short int *lowcl, short int *highcl,
99 float *outval, float *outvalc, float *outvalv);
100static cpl_table *mknewtab(casu_tfits *cat, dstrct *ddcur, dstrct *ddref);
101static cpl_array *transinit(void);
102static void fileptrs_close(int nim, dstrct *fileptrs);
103static void tidy(int nim, dstrct *fileptrs, cpl_wcs *outwcs);
104
105
108/*---------------------------------------------------------------------------*/
185/*---------------------------------------------------------------------------*/
186
187extern int casu_imstack(casu_fits **inf, casu_fits **inconf, casu_fits **invar,
188 casu_tfits **cats, int nimages, int nconfs, float lthr,
189 float hthr, int method, int seeing, int fast, int unmap,
190 const char *expkey, casu_fits **out, casu_fits **outc,
191 casu_fits **outv, int *status) {
192 int i,nm;
193 dstrct *dd;
194 float expref,exposure,xoff,yoff,sumweights;
195 casu_fits *outfits=NULL,*outconf=NULL,*outvar=NULL;
196 casu_tfits *catref;
197 cpl_table *newtab,*mtab;
198 cpl_propertylist *phu;
199 cpl_wcs *outwcs = NULL;
200 dstrct *fileptrs=NULL;
201 const char *fctid = "casu_imstack";
202
203 /* Inherited status */
204
205 *out = NULL;
206 *outc = NULL;
207 *outv = NULL;
208 if (*status != CASU_OK)
209 return(*status);
210
211 /* Is there any point in being here? */
212
213 if (nimages == 0) {
214 cpl_msg_error(fctid,"No input files to combine");
215 tidy(nimages,fileptrs,outwcs);
216 FATAL_ERROR
217 } else if (nimages == 1) {
218 *out = casu_fits_duplicate(inf[0]);
219 if (invar != NULL)
220 *outv = casu_fits_duplicate(invar[0]);
221 *outc = casu_fits_duplicate(inconf[0]);
222 casu_prov(casu_fits_get_ehu(*out),inf,nimages,1);
223 return(*status);
224 }
225
226 /* Allocate file struct array and fill in some values */
227
228 fileptrs = cpl_malloc(nimages*sizeof(dstrct));
229 if (cpl_propertylist_has(casu_fits_get_phu(inf[0]),expkey)) {
230 exposure = (float)cpl_propertylist_get_double(casu_fits_get_phu(inf[0]),
231 expkey);
232 } else {
233 exposure = 0.5;
234 }
235 expref = max(0.5,exposure);
236 catref = NULL;
237 for (i = 0; i < nimages; i++) {
238 dd = fileptrs + i;
239 dd->fname = inf[i];
240
241 /* Get the confidence map and it's header scale */
242
243 dd->confcopy = 0;
244 if (nconfs == 0) {
245 dd->conf = NULL;
246 } else if (nconfs == 1) {
247 dd->conf = inconf[0];
248 dd->confcopy = 1;
249 } else {
250 dd->conf = inconf[i];
251 }
252
253 /* Load up the variance arrays if they exist */
254
255 if (invar != NULL)
256 dd->var = invar[i];
257 else
258 dd->var = NULL;
259
260 /* Get the image size */
261
262 dd->nx = (int)cpl_image_get_size_x(casu_fits_get_image(dd->fname));
263 dd->ny = (int)cpl_image_get_size_y(casu_fits_get_image(dd->fname));
264 if (unmap)
265 casu_fits_unload_im(dd->fname);
266
267 /* Get the RA/Dec range covered by this image (only needed for
268 slow algorithm) */
269
270 if (fast == 0) {
271 *status = CASU_OK;
272 (void)casu_coverage(casu_fits_get_ehu(dd->fname),0,&(dd->ramin),
273 &(dd->ramax),&(dd->decmin),&(dd->decmax),
274 status);
275 } else {
276 dd->ramin = 0.0;
277 dd->ramax = 0.0;
278 dd->decmin = 0.0;
279 dd->decmax = 0.0;
280 }
281
282 /* Get exposure time info */
283
284 phu = casu_fits_get_phu(dd->fname);
285 if (cpl_propertylist_has(phu,expkey)) {
286 exposure = (float)cpl_propertylist_get_double(phu,expkey);
287 exposure = max(0.5,exposure);
288 } else {
289 exposure = 0.5;
290 }
291 dd->expscale = exposure/expref;
292
293 /* Zero the bad pixel mask */
294
295 dd->bpm = cpl_calloc((dd->nx)*(dd->ny),sizeof(unsigned char));
296
297 /* Read the WCS from the header */
298
299 dd->vwcs = cpl_wcs_new_from_propertylist(casu_fits_get_ehu(inf[i]));
300
301 /* Double check to make sure the confidence maps and images have the
302 same dimensions */
303
304 if ((int)cpl_image_get_size_x(casu_fits_get_image(dd->conf)) != dd->nx ||
305 (int)cpl_image_get_size_y(casu_fits_get_image(dd->conf)) != dd->ny) {
306 cpl_msg_error(fctid,"Image %s and Confidence map %s don't match",
307 casu_fits_get_fullname(dd->fname),
308 casu_fits_get_fullname(dd->conf));
309 tidy(nimages,fileptrs,outwcs);
310 FATAL_ERROR
311 }
312
313 /* Double check to make sure the variance maps and images have the
314 same dimensions */
315
316 if (dd->var != NULL) {
317 if ((int)cpl_image_get_size_x(casu_fits_get_image(dd->var)) != dd->nx ||
318 (int)cpl_image_get_size_y(casu_fits_get_image(dd->var)) != dd->ny) {
319 cpl_msg_error(fctid,"Image %s and Variance map %s don't match",
320 casu_fits_get_fullname(dd->fname),
321 casu_fits_get_fullname(dd->var));
322 tidy(nimages,fileptrs,outwcs);
323 FATAL_ERROR
324 }
325 if (unmap)
326 casu_fits_unload_im(dd->var);
327 }
328
329 /* Get a rough shift between the frames and the first one */
330
331 diffxy(fileptrs->vwcs,dd->vwcs,&(dd->xoff),&(dd->yoff));
332
333 /* Do we have catalogues? If so, then work out the transformation
334 corrections from them */
335
336 if (cats != NULL) {
337 if (cats[i] != NULL)
338 catref = cats[i];
339 if (i == 0 || cats[i] == NULL) {
340 dd->trans = transinit();
341 dd->invtrans = transinit();
342 } else {
343 newtab = mknewtab(cats[i],dd,fileptrs);
344 (void)casu_matchxy(newtab,casu_tfits_get_table(catref),
345 1024,&xoff,&yoff,&nm,&mtab,status);
346 if (*status != CASU_OK) {
347 dd->trans = transinit();
348 *status = CASU_OK;
349 } else {
350 (void)casu_platexy(mtab,6,&(dd->trans),status);
351 freetable(mtab);
352 (void)casu_matchxy(casu_tfits_get_table(catref),newtab,
353 1024,&xoff,&yoff,&nm,&mtab,status);
354 (void)casu_platexy(mtab,6,&(dd->invtrans),status);
355 }
356 freetable(mtab);
357 freetable(newtab);
358 }
359 } else {
360 dd->trans = NULL;
361 dd->invtrans = NULL;
362 }
363 }
364
365 /* Get the background levels in the overlap regions */
366
367 backgrnd_ov(nimages,fileptrs,unmap);
368
369 /* Do the seeing weighting if you want to */
370
371 if (seeing)
372 seeing_wt(nimages,fileptrs);
373
374 /* Get a sum of weights */
375
376 sumweights = 0.0;
377 for (i = 0; i < nimages; i++)
378 sumweights += (fileptrs+i)->weight;
379
380 /* Create the output file */
381
382 output_files(nimages,fileptrs,&outfits,&outconf,&outvar,&outwcs);
383
384 /* Do the stacking for the fast routines */
385
386 if (fast == 1) {
387
388 /* Our first job is to do a nearest neighbour stack in order
389 to find the objects that will be rejected. If we only want
390 the NN algorithm, then we're done. Otherwise we can use the
391 rejection information to restack using any other algorithm */
392
393 stack_nn(nimages,fileptrs,lthr,hthr,sumweights,outfits,outconf,
394 outvar,outwcs);
395
396 /* Switch for final stacking method */
397
398 switch (method) {
399 case 0:
400 break;
401 case 1:
402 stack_lin(nimages,fileptrs,sumweights,outfits,outconf,outvar,outwcs);
403 break;
404 default:
405 break;
406 }
407
408 /* Otherwise use the slow routine */
409
410 } else {
411 stackslow(nimages,fileptrs,lthr,hthr,sumweights,outfits,outconf,
412 outvar,outwcs,method);
413 }
414
415 /* Add the provenance cards to the header */
416
417 casu_prov(casu_fits_get_ehu(outfits),inf,nimages,1);
418
419 /* Assign some output info, tidy and get out of here */
420
421 *out = outfits;
422 *outc = outconf;
423 *outv = outvar;
424 tidy(nimages,fileptrs,outwcs);
425 return(CASU_OK);
426}
427
428/*---------------------------------------------------------------------------*/
465/*---------------------------------------------------------------------------*/
466
467static void stackslow(int nimages, dstrct *fileptrs, float lsig, float hsig,
468 float sumweights, casu_fits *outfits, casu_fits *outconf,
469 casu_fits *outvar, cpl_wcs *outwcs, int dolin) {
470 float scale_avvar,*outdata,*outdatav,lclip,hclip,*data;
471 float *vdata = NULL; /* return value for cpl_image_get_data_float on error */
472 float *workbuf,*workbufv,*workbufw,*workbufc,*dbuf,*wbuf,*cbuf;
473 float outval,outvalc,outvalv,avlev,avvar,dx,dy,value,cval,*confwork;
474 float med,wt = 0.,w[4],renorm,*vbuf,vval,*workdatav,*workdatav2,*workdatav3;
475 double yout,xin,yin,oxf,oyf,lra1,lra2,ldec1,ldec2;
476 int *nbuf,*iloc,jout,bufloc_row,bufloc_before,iim,i,j,*outdatac;
477 int gotit,iout,ixin,iyin,ind,ind1,ox,oy,bufloc,ncontrib,*lbuf,nin;
478 int xmin,xmax,ymin,ymax,nz,nx,*bloc,*bbuf,m,nn,ii,jj,*cdata;
479 short int clipped_low,clipped_high;
480 long nxo,nyo,npts;
481 cpl_image *im,*cim,*vim;
482 unsigned char *clipmap,*id,*ibuf,**bpms,*bbpm;
483 dstrct *dd;
484 cpl_matrix *xyin,*xyin2,*xyout;
485
486 outdatav = NULL; /* Remove falsepositive compiler warning */
487 workdatav = NULL; /* Remove falsepositive compiler warning from freespace macro*/
488 workdatav2 = NULL; /* Remove falsepositive compiler warning from freespace macro*/
489 workdatav3 = NULL; /* Remove falsepositive compiler warning from freespace macro*/
490
491 /* Useful constant */
492
493 scale_avvar = sqrt(CPL_MATH_PI/2.0)/9.0;
494
495 /* Output file info */
496
497 outdata = cpl_image_get_data_float(casu_fits_get_image(outfits));
498 if (outvar != NULL)
499 outdatav = cpl_image_get_data_float(casu_fits_get_image(outvar));
500 outdatac = cpl_image_get_data_int(casu_fits_get_image(outconf));
501 nxo = cpl_image_get_size_x(casu_fits_get_image(outfits));
502 nyo = cpl_image_get_size_y(casu_fits_get_image(outfits));
503
504 /* Rejection thresholds */
505
506 lclip = fileptrs->sky - lsig*(fileptrs->noise);
507 hclip = fileptrs->sky + hsig*(fileptrs->noise);
508
509 /* Get some workspace and initialise it */
510
511 nz = nimages + 2;
512 workbuf = cpl_malloc(2*nxo*nz*sizeof(float));
513 if (outvar != NULL)
514 workbufv = cpl_malloc(2*nxo*nz*sizeof(float));
515 else
516 workbufv = NULL;
517 workbufw = cpl_malloc(2*nxo*nz*sizeof(float));
518 workbufc = cpl_malloc(2*nxo*nz*sizeof(float));
519 confwork = cpl_calloc(nxo*nyo,sizeof(float));
520 id = cpl_malloc(2*nxo*nz*sizeof(unsigned char));
521 nbuf = cpl_malloc(2*nxo*sizeof(int));
522 iloc = cpl_malloc(2*nxo*nz*sizeof(int));
523 bloc = cpl_malloc(2*nxo*nz*sizeof(int));
524 clipmap = cpl_calloc(2*nxo*nz,sizeof(unsigned char));
525 bpms = cpl_calloc(nimages,sizeof(unsigned char *));
526 dbuf = cpl_malloc(2*nz*sizeof(float));
527 wbuf = cpl_malloc(2*nz*sizeof(float));
528 if (outvar != NULL)
529 vbuf = cpl_malloc(2*nz*sizeof(float));
530 else
531 vbuf = NULL;
532 cbuf = cpl_malloc(2*nz*sizeof(float));
533 ibuf = cpl_malloc(2*nz*sizeof(unsigned char));
534 lbuf = cpl_malloc(2*nz*sizeof(int));
535 bbuf = cpl_malloc(2*nz*sizeof(int));
536
537 /* Get a bad pixel mask for each input file */
538
539 for (i = 0; i < nimages; i++)
540 bpms[i] = (fileptrs+i)->bpm;
541
542 /* Loop for each output row */
543
544 for (jout = 0; jout < nyo; jout++) {
545 bufloc_row = (jout % 2 ? 0 : 1);
546 bufloc_before = (! bufloc_row);
547 yout = (double)(jout+1);
548
549 /* Zero the counter for this row */
550
551 (void)memset(nbuf+bufloc_row*nxo,0,nxo*sizeof(int));
552
553 /* Create an input matrix of positions and work out the
554 coverage of the output line */
555
556 nin = 2*nxo;
557 xyin = cpl_matrix_new((cpl_size)nin,2);
558 nin = 0;
559 for (i = -1; i <= 1; i+= 2) {
560 for (iout = 0; iout < nxo; iout++) {
561 cpl_matrix_set(xyin,nin,0,(double)(iout+1));
562 cpl_matrix_set(xyin,nin,1,yout+0.5*(double)i);
563 nin++;
564 }
565 }
566 linecoverage(outwcs,xyin,&lra1,&lra2,&ldec1,&ldec2);
567
568 /* Loop for each image and find the range of input pixels for
569 this image that will affect this output row */
570
571 for (iim = 0; iim < nimages; iim++) {
572
573 /* Initialise a few things */
574
575 dd = fileptrs + iim;
576 wt = dd->weight;
577 gotit = 0;
578
579 /* If this line has nothing in common with the current
580 image, then move to the next image */
581
582 if (ldec2 < dd->decmin || ldec1 > dd->decmax ||
583 lra2 < dd->ramin || lra1 > dd->ramax)
584 continue;
585
586 /* Get the position matrix */
587
588 inloc(outwcs,xyin,dd->vwcs,dd->trans,&xyout);
589
590 /* Get the xy limits now */
591
592 xmin = 1000000000;
593 xmax = -1000000000;
594 ymin = 1000000000;
595 ymax = -1000000000;
596 nin = 0;
597 for (i = -1; i <= 1; i += 2) {
598 for (iout = 0; iout < nxo; iout++) {
599 xin = cpl_matrix_get(xyout,nin,0);
600 yin = cpl_matrix_get(xyout,nin,1);
601 nin++;
602 ixin = casu_nint(xin);
603 iyin = casu_nint(yin);
604 if (ixin < 1 || iyin < 1 || ixin > dd->nx ||
605 iyin > dd->ny)
606 continue;
607 xmin = min(xmin,ixin);
608 xmax = max(xmax,ixin);
609 ymin = min(ymin,iyin);
610 ymax = max(ymax,iyin);
611 gotit = 1;
612 }
613 }
614 cpl_matrix_delete(xyout);
615
616 /* If this makes a contribution, then get the data subset */
617
618 if (! gotit)
619 continue;
620 im = cpl_image_load_window(casu_fits_get_filename(dd->fname),
621 CPL_TYPE_FLOAT,0,
622 (cpl_size)casu_fits_get_nexten(dd->fname),
623 (cpl_size)xmin,(cpl_size)ymin,
624 (cpl_size)xmax,(cpl_size)ymax);
625 nx = cpl_image_get_size_x(im);
626 data = cpl_image_get_data_float(im);
627 cim = cpl_image_load_window(casu_fits_get_filename(dd->conf),
628 CPL_TYPE_INT,0,
629 (cpl_size)casu_fits_get_nexten(dd->fname),
630 (cpl_size)xmin,(cpl_size)ymin,
631 (cpl_size)xmax,(cpl_size)ymax);
632 cdata = cpl_image_get_data_int(cim);
633 if (outvar != NULL) {
634 vim = cpl_image_load_window(casu_fits_get_filename(dd->var),
635 CPL_TYPE_FLOAT,0,
636 (cpl_size)casu_fits_get_nexten(dd->fname),
637 (cpl_size)xmin,(cpl_size)ymin,
638 (cpl_size)xmax,(cpl_size)ymax);
639 vdata = cpl_image_get_data_float(vim);
640 }
641
642 /* Get the position matrix */
643
644 nin = (ymax - ymin + 1)*(xmax - xmin + 1);
645 xyin2 = cpl_matrix_new((cpl_size)nin,2);
646 nin = 0;
647 for (j = ymin; j <= ymax; j++) {
648 for (i = xmin; i <= xmax; i++) {
649 cpl_matrix_set(xyin2,nin,0,(double)i);
650 cpl_matrix_set(xyin2,nin,1,(double)j);
651 nin++;
652 }
653 }
654 outloc(dd->vwcs,xyin2,outwcs,dd->trans,&xyout);
655 cpl_matrix_delete(xyin2);
656
657 /* Ok, distribute the data into the buffer */
658
659 nin = 0;
660 for (j = ymin; j <= ymax; j++) {
661 ind = (j - ymin)*nx;
662 for (i = xmin; i <= xmax; i++) {
663 ind1 = ind + i - xmin;
664 oxf = cpl_matrix_get(xyout,nin,0);
665 oyf = cpl_matrix_get(xyout,nin,1);
666 nin++;
667 ox = casu_nint(oxf) - 1;
668 oy = casu_nint(oyf) - 1;
669 if (ox < 0 || ox >= nxo || oy < 0 || oy >= nyo ||
670 cdata[ind1] == 0 || oy != jout)
671 continue;
672 bufloc = (bufloc_row*nxo + ox)*nz +
673 nbuf[ox+bufloc_row*nxo];
674 value = data[ind1]/(dd->expscale) + dd->skydiff;
675 workbuf[bufloc] = value;
676 if (outvar != NULL)
677 workbufv[bufloc] = vdata[ind1]/powf(dd->expscale,2.0);
678 workbufw[bufloc] = wt;
679 workbufc[bufloc] = (float)cdata[ind1];
680 id[bufloc] = iim;
681 nbuf[ox+bufloc_row*nxo] += 1;
682 iloc[bufloc] = ind1;
683 bloc[bufloc] = (j-1)*(dd->nx) + i - 1;
684 bpms[iim][bloc[bufloc]] = 0;
685 }
686 }
687 cpl_matrix_delete(xyout);
688
689 /* Do some intermediate tidying */
690
691 freeimage(im);
692 freeimage(cim);
693 if (outvar != NULL)
694 freeimage(vim);
695 } /* End of image loop */
696 cpl_matrix_delete(xyin);
697
698 /* Loop through the work buffer and find an initial estimate
699 for the current output row */
700
701 for (iout = 0; iout < nxo; iout++) {
702 clipmap[iout+bufloc_row*nxo] = 0;
703 bufloc = (bufloc_row*nxo + iout)*nz;
704 ncontrib = nbuf[iout+bufloc_row*nxo];
705 if (ncontrib == 0) {
706 outdata[jout*nxo + iout] = fileptrs->sky;
707 outdatac[jout*nxo + iout] = 0;
708 continue;
709 }
710
711 /* Put some stuff in the buffers for the averaging routine */
712
713 for (i = 0; i < ncontrib; i++) {
714 dbuf[i] = workbuf[bufloc+i];
715 wbuf[i] = workbufw[bufloc+i];
716 if (outvar != NULL)
717 vbuf[i] = workbufv[bufloc+i];
718 cbuf[i] = workbufc[bufloc+i];
719 ibuf[i] = id[bufloc+i];
720 lbuf[i] = iloc[bufloc+i];
721 bbuf[i] = bloc[bufloc+i];
722 }
723
724 /* Do the averages with nominal clipping */
725
726 do_averages(ncontrib,dbuf,vbuf,wbuf,cbuf,ibuf,lclip,hclip,hsig,0.0,
727 fileptrs,sumweights,&clipped_low,&clipped_high,&outval,
728 &outvalc,&outvalv);
729
730 /* Store these away */
731
732 outdata[jout*nxo + iout] = outval;
733 confwork[jout*nxo + iout] = outvalc;
734 if (outvar != NULL)
735 outdatav[jout*nxo + iout] = outvalv;
736 clipmap[bufloc_row*nxo + iout] = (clipped_high >= 0);
737 if (clipped_low >= 0)
738 bpms[ibuf[clipped_low]][bbuf[clipped_low]] = 2;
739 if (clipped_high >= 0)
740 bpms[ibuf[clipped_high]][bbuf[clipped_high]] = 1;
741 }
742
743 /* Look at the non-edge pixels and see if the local variation can
744 show whether any clipping that has been done is justified. NB:
745 we are looking at the _previous_ row here */
746
747 if (jout >= 2) {
748 for (iout = 1; iout < nxo-1; iout++) {
749 if (! clipmap[bufloc_before*nxo + iout])
750 continue;
751 ncontrib = nbuf[bufloc_before*nxo+iout];
752 if (ncontrib == 0)
753 continue;
754 ind = (jout-1)*nxo + iout;
755 bufloc = (bufloc_before*nxo + iout)*nz;
756 avlev = 0.0;
757 for (j = -1; j <= 1; j++)
758 for (i = -1; i <= 1; i++)
759 avlev += outdata[ind + j*nxo + i];
760 avlev /= 9.0;
761 avvar = 0.0;
762 for (j = -1; j <= 1; j++)
763 for (i = -1; i <= 1; i++)
764 avvar += (float)fabs(avlev - outdata[ind + j*nxo + i]);
765 avvar *= scale_avvar;
766 if (avlev <= hclip && avvar <= 2.0*(fileptrs->noise))
767 continue;
768
769 /* Put some stuff in the buffers for the averaging routine */
770
771 for (i = 0; i < ncontrib; i++) {
772 dbuf[i] = workbuf[bufloc+i];
773 if (outvar != NULL)
774 vbuf[i] = workbufv[bufloc+i];
775 wbuf[i] = workbufw[bufloc+i];
776 cbuf[i] = workbufc[bufloc+i];
777 ibuf[i] = id[bufloc+i];
778 lbuf[i] = iloc[bufloc+i];
779 bbuf[i] = bloc[bufloc+i];
780 bpms[ibuf[i]][bbuf[i]] = 0;
781 }
782
783 /* Redo the averages with the nominal clipping */
784
785 do_averages(ncontrib,dbuf,vbuf,wbuf,cbuf,ibuf,lclip,hclip,hsig,
786 avvar,fileptrs,sumweights,&clipped_low,
787 &clipped_high,&outval,&outvalc,&outvalv);
788
789 /* Store these away */
790
791 outdata[ind] = outval;
792 if (outvar != NULL)
793 outdatav[ind] = outvalv;
794 confwork[ind] = outvalc;
795 if (clipped_low >= 0)
796 bpms[ibuf[clipped_low]][bbuf[clipped_low]] = 2;
797 if (clipped_high >= 0)
798 bpms[ibuf[clipped_high]][bbuf[clipped_high]] = 1;
799 }
800 }
801 }
802
803 /* Do we want to do the bi-linear interpolation now? If so,
804 the first thing to do is to throw away the NN result but
805 keep the information about which input pixels were rejected */
806
807 if (dolin) {
808 if (outvar != NULL) {
809 workdatav = cpl_malloc(nxo*nyo*sizeof(float));
810 workdatav2 = cpl_malloc(nxo*nyo*sizeof(float));
811 workdatav3 = cpl_calloc(nxo*nyo,sizeof(float));
812 }
813 memset(outdata,0,nxo*nyo*sizeof(float));
814 memset(confwork,0,nxo*nyo*sizeof(float));
815 if (outvar != NULL)
816 memset(outdatav,0,nxo*nyo*sizeof(float));
817
818 /* Get the data for the current input image */
819
820 npts = nxo*nyo;
821 for (iim = 0; iim < nimages; iim++) {
822 if (outvar != NULL) {
823 memset(workdatav,0,nxo*nyo*sizeof(float));
824 memset(workdatav2,0,nxo*nyo*sizeof(float));
825 }
826 dd = fileptrs + iim;
827 im = cpl_image_load(casu_fits_get_filename(dd->fname),
828 CPL_TYPE_FLOAT,0,
829 (cpl_size)casu_fits_get_nexten(dd->fname));
830 data = cpl_image_get_data_float(im);
831 cim = cpl_image_load(casu_fits_get_filename(dd->conf),
832 CPL_TYPE_INT,0,
833 (cpl_size)casu_fits_get_nexten(dd->fname));
834 cdata = cpl_image_get_data_int(cim);
835 if (outvar != NULL) {
836 vim = cpl_image_load(casu_fits_get_filename(dd->var),
837 CPL_TYPE_FLOAT,0,
838 (cpl_size)casu_fits_get_nexten(dd->fname));
839 vdata = cpl_image_get_data_float(vim);
840 }
841
842 /* Get position matrix */
843
844 nin = (dd->ny)*(dd->nx);
845 xyin = cpl_matrix_new((cpl_size)nin,2);
846 nin = 0;
847 for (j = 1; j <= dd->ny; j++) {
848 for (i = 1; i <= dd->nx; i++) {
849 cpl_matrix_set(xyin,nin,0,(double)i);
850 cpl_matrix_set(xyin,nin,1,(double)j);
851 nin++;
852 }
853 }
854 outloc(dd->vwcs,xyin,outwcs,dd->trans,&xyout);
855 cpl_matrix_delete(xyin);
856
857 /* Now distribute the data */
858
859 nin = -1;
860 for (j = 1; j <= dd->ny; j++) {
861 ind = (j - 1)*(dd->nx);
862 for (i = 1; i <= dd->nx; i++) {
863 nin++;
864 ind1 = ind + i - 1;
865 if (bpms[iim][ind1] != 0 || cdata[ind1] == 0)
866 continue;
867 oxf = cpl_matrix_get(xyout,nin,0);
868 oyf = cpl_matrix_get(xyout,nin,1);
869 if (oxf < 0.0 || oyf < 0.0)
870 continue;
871 ox = (int)oxf - 1;
872 oy = (int)oyf - 1;
873 if (ox < 0 || ox >= nxo || oy < 0 || oy >= nyo)
874 continue;
875 dx = oxf - (float)(ox + 1);
876 dy = oyf - (float)(oy + 1);
877 w[0] = (1.0-dx)*(1.0-dy);
878 w[1] = dx*(1.0-dy);
879 w[2] = dy*(1.0-dx);
880 w[3] = dx*dy;
881 m = 0;
882 value = data[ind1]/(dd->expscale) + dd->skydiff;
883 if (outvar != NULL)
884 vval = vdata[ind1]/powf(dd->expscale,2.0);
885 cval = (float)cdata[ind1];
886 for (jj = 0; jj <= 1; jj++) {
887 if (oy+jj >= nyo)
888 continue;
889 for (ii = 0; ii <= 1; ii++) {
890 nn = (oy+jj)*nxo + ox + ii;
891 if (nn >= npts || ox+ii >= nxo)
892 continue;
893 outdata[nn] += w[m]*value*wt*cval;
894 confwork[nn] += w[m]*wt*cval;
895 if (outvar != NULL) {
896 workdatav[nn] += vval*w[m]*cval;
897 workdatav2[nn] += w[m]*cval;
898 }
899 m++;
900 }
901 }
902 }
903 }
904 cpl_matrix_delete(xyout);
905 freeimage(im);
906 freeimage(cim);
907 if (outvar != NULL) {
908 freeimage(vim);
909 for (ii = 0; ii < npts; ii++) {
910 if (workdatav2[ii] != 0.0) {
911 outdatav[ii] += wt*wt*workdatav[ii]/workdatav2[ii];
912 workdatav3[ii] += wt;
913 }
914 }
915 }
916 }
917
918 /* Now average the output data */
919
920 for (jout = 0; jout < nyo; jout++) {
921 for (iout = 0; iout < nxo; iout++) {
922 ind = jout*nxo + iout;
923 if (confwork[ind] != 0.0) {
924 outdata[ind] /= confwork[ind];
925 if (outvar != NULL) {
926 if (workdatav3[ind] != 0.0)
927 outdatav[ind] /= powf(workdatav3[ind],2.0);
928 else
929 outdatav[ind] = 0.0;
930 }
931 } else {
932 outdata[ind] = fileptrs->sky;
933 if (outvar != NULL)
934 outdatav[ind] = 0.0;
935 }
936 confwork[ind] /= sumweights;
937 }
938 }
939 if (outvar != NULL) {
940 freespace(workdatav);
941 freespace(workdatav2);
942 freespace(workdatav3);
943 }
944 }
945
946 /* Normalise the confidence */
947
948 npts = nxo*nyo;
949 bbpm = cpl_calloc(npts,sizeof(*bbpm));
950 for (i = 0; i < npts; i++)
951 bbpm[i] = (confwork[i] == 0.0);
952 med = casu_med(confwork,bbpm,npts);
953 freespace(bbpm);
954 renorm = 100.0/med;
955 for (i = 0; i < npts; i++) {
956 confwork[i] *= renorm;
957 outdatac[i] = max(0,min(1000,(int)(confwork[i]+0.5)));
958 }
959
960 /* Tidy up */
961
962 freespace(bpms);
963 freespace(workbuf);
964 freespace(workbufw);
965 freespace(workbufv);
966 freespace(workbufc);
967 freespace(confwork);
968 freespace(id);
969 freespace(nbuf);
970 freespace(iloc);
971 freespace(bloc);
972 freespace(clipmap);
973 freespace(dbuf);
974 freespace(vbuf);
975 freespace(wbuf);
976 freespace(cbuf);
977 freespace(ibuf);
978 freespace(lbuf);
979 freespace(bbuf);
980}
981
982/*---------------------------------------------------------------------------*/
1009/*---------------------------------------------------------------------------*/
1010
1011static void linecoverage(cpl_wcs *outwcs, cpl_matrix *xyin, double *lra1,
1012 double *lra2, double *ldec1, double *ldec2) {
1013 int i,first_quad,fourth_quad,nr;
1014 cpl_matrix *radec;
1015 cpl_array *status;
1016 double cur_ra,cur_dec,*data,max_1q,min_4q;
1017
1018 /* Do the physical -> world conversion */
1019
1020 cpl_wcs_convert(outwcs,xyin,&radec,&status,CPL_WCS_PHYS2WORLD);
1021 cpl_array_delete(status);
1022
1023 /* Work out the min and max of the RA and Dec */
1024
1025 *lra1 = 370.0;
1026 *lra2 = -370.0;
1027 *ldec1 = 95;
1028 *ldec2 = -95.0;
1029 first_quad = 0;
1030 fourth_quad = 0;
1031 data = cpl_matrix_get_data(radec);
1032 nr = cpl_matrix_get_nrow(radec);
1033 min_4q = -91.0;
1034 max_1q = data[0]; /* Added as max_1q needs to be initialized */
1035 for (i = 0; i < nr; i++) {
1036 cur_ra = data[2*i];
1037 cur_dec = data[2*i+1];
1038 if (cur_ra >= 0.0 && cur_ra <= 90.0) {
1039 first_quad = 1;
1040 max_1q = max(cur_ra,max_1q);
1041 } else if (cur_ra >= 270.0 && cur_ra < 360.0) {
1042 fourth_quad = 1;
1043 min_4q = min((cur_ra-360.0),min_4q);
1044 }
1045 *lra1 = min(*lra1,cur_ra);
1046 *lra2 = max(*lra2,cur_ra);
1047 *ldec1 = min(*ldec1,cur_dec);
1048 *ldec2 = max(*ldec2,cur_dec);
1049 }
1050 cpl_matrix_delete(radec);
1051
1052 /* Now look to see if you had RA values in both the first and
1053 fourth quadrants. If you have, then make the minimum RA a negative
1054 value. This will be the signal to the caller that you have a
1055 wraparound... */
1056
1057 if (first_quad && fourth_quad) {
1058 *lra1 = min_4q;
1059 *lra2 = max_1q;
1060 }
1061}
1062
1063/*---------------------------------------------------------------------------*/
1097/*---------------------------------------------------------------------------*/
1098
1099static int stack_nn(int nim, dstrct *fileptrs, float lsig, float hsig,
1100 float sumweights, casu_fits *outfits, casu_fits *outconf,
1101 casu_fits *outvar, cpl_wcs *outwcs) {
1102 float *workbuf,*workbufc,*workbufw,lclip,hclip,*data,wt,med,outvalc;
1103 float value,*dbuf,*wbuf,*cbuf,outval,avlev,avvar,*confwork,*outdata;
1104 float scale_avvar,renorm,*vbuf,*outdatav,*workbufv,*vdata=NULL,outvalv;
1105 double oxf,oyf;
1106 long *iloc,i,npts,ix,iy,ind1,ind,ox,oy,nn,bufloc,ncontrib,*lbuf;
1107 long j,nx,ny,nz,nxo,nyo,nall;
1108 int *outdatac,*cdata,nin;
1109 short int clipped_low,clipped_high;
1110 dstrct *dd;
1111 unsigned char *id,*clipmap,*cc,*cclast,*ibuf,*bbpm;
1112 short int *nbuf;
1113 cpl_image *outim,*outimc,*outimv;
1114 cpl_matrix *xyin,*xyout;
1115
1116 outdatav = NULL; /* Remove falsepositive compiler warning */
1117
1118 /* Useful constant */
1119
1120 scale_avvar = sqrt(CPL_MATH_PI/2.0)/9.0;
1121
1122 /* Output file info */
1123
1124 outim = casu_fits_get_image(outfits);
1125 outimc = casu_fits_get_image(outconf);
1126 if (outvar != NULL)
1127 outimv = casu_fits_get_image(outvar);
1128 nxo = (int)cpl_image_get_size_x(outim);
1129 nyo = (int)cpl_image_get_size_y(outim);
1130
1131 /* Get workspace for collecting info over the whole output map */
1132
1133 nz = nim + 2;
1134 nall = (long)nz*(long)nxo*(long)nyo;
1135 workbuf = cpl_calloc(nall,sizeof(*workbuf));
1136 if (outvar != NULL)
1137 workbufv = cpl_calloc(nall,sizeof(*workbufv));
1138 else
1139 workbufv = NULL;
1140 workbufc = cpl_calloc(nall,sizeof(*workbufc));
1141 workbufw = cpl_calloc(nall,sizeof(*workbufw));
1142 id = cpl_calloc(nall,sizeof(*id));
1143 iloc = cpl_calloc(nall,sizeof(*iloc));
1144 nbuf = cpl_calloc(nxo*nyo,sizeof(*nbuf));
1145 confwork = cpl_calloc(nxo*nyo,sizeof(*confwork));
1146
1147 /* Workspace for marking output rows with clipped pixels */
1148
1149 clipmap = cpl_calloc(2*nxo,sizeof(*clipmap));
1150
1151 /* Buffers for averaging an individual pixel */
1152
1153 dbuf = cpl_calloc(2*nz,sizeof(*dbuf));
1154 if (outvar != NULL)
1155 vbuf = cpl_calloc(2*nz,sizeof(*vbuf));
1156 else
1157 vbuf = NULL;
1158 wbuf = cpl_calloc(2*nz,sizeof(*wbuf));
1159 cbuf = cpl_calloc(2*nz,sizeof(*cbuf));
1160 ibuf = cpl_calloc(2*nz,sizeof(*ibuf));
1161 lbuf = cpl_calloc(2*nz,sizeof(*lbuf));
1162
1163 /* Get output buffers */
1164
1165 outdata = cpl_image_get_data_float(outim);
1166 outdatac = cpl_image_get_data_int(outimc);
1167 if (outvar != NULL)
1168 outdatav = cpl_image_get_data_float(outimv);
1169
1170 /* Rejection thresholds */
1171
1172 lclip = fileptrs->sky - lsig*(fileptrs->noise);
1173 hclip = fileptrs->sky + hsig*(fileptrs->noise);
1174
1175 /* Ok, loop for each of the input files and collect all the information
1176 for each pixel */
1177
1178 for (i = 0; i < nim; i++) {
1179 dd = fileptrs + i;
1180 nx = dd->nx;
1181 ny = dd->ny;
1182 npts = nx*ny;
1183 data = cpl_image_get_data_float(casu_fits_get_image(dd->fname));
1184 cdata = cpl_image_get_data_int(casu_fits_get_image(dd->conf));
1185 if (outvar != NULL)
1186 vdata = cpl_image_get_data_float(casu_fits_get_image(dd->var));
1187 wt = dd->weight;
1188
1189 /* Get the position matrix */
1190
1191 nin = nx*ny;
1192 xyin = cpl_matrix_new((cpl_size)nin,2);
1193 nin = 0;
1194 for (iy = 1; iy <= ny; iy++) {
1195 for (ix = 1; ix <= nx; ix++) {
1196 cpl_matrix_set(xyin,nin,0,(double)ix);
1197 cpl_matrix_set(xyin,nin,1,(double)iy);
1198 nin++;
1199 }
1200 }
1201 outloc(dd->vwcs,xyin,outwcs,dd->trans,&xyout);
1202 cpl_matrix_delete(xyin);
1203
1204 /* Loop for each pixel in the data array and work out where it belongs
1205 in the output array */
1206
1207 nin = 0;
1208 for (iy = 0; iy < ny; iy++) {
1209 ind1 = iy*nx;
1210 for (ix = 0; ix < nx; ix++) {
1211 ind = ind1 + ix;
1212 oxf = cpl_matrix_get(xyout,nin,0);
1213 oyf = cpl_matrix_get(xyout,nin,1);
1214 nin++;
1215 ox = casu_nint(oxf) - 1;
1216 oy = casu_nint(oyf) - 1;
1217 if (ox < 0 || ox >= nxo || oy < 0 || oy >= nyo ||
1218 cdata[ind] == 0)
1219 continue;
1220 nn = oy*nxo + ox;
1221 bufloc = oy*nxo*nz + ox*nz + nbuf[nn];
1222 value = data[ind]/(dd->expscale) + dd->skydiff;
1223 workbuf[bufloc] = value;
1224 if (outvar != NULL && vdata != NULL)
1225 workbufv[bufloc] = vdata[ind]/((dd->expscale)*(dd->expscale));
1226 workbufw[bufloc] = wt;
1227 workbufc[bufloc] = (float)cdata[ind];
1228 id[bufloc] = i;
1229 nbuf[nn] += 1;
1230 iloc[bufloc] = ind;
1231 }
1232 }
1233 cpl_matrix_delete(xyout);
1234 }
1235
1236 /* Now loop through the output arrays and form an initial estimate */
1237
1238 for (iy = 0; iy < nyo; iy++) {
1239 ind1 = iy*nxo;
1240 cc = clipmap + (iy % 2)*nxo;
1241 cclast = clipmap + ((iy-1) % 2)*nxo;
1242 for (ix = 0; ix < nxo; ix++) {
1243 ind = ind1 + ix;
1244 bufloc = iy*nxo*nz + ix*nz;
1245 ncontrib = nbuf[ind];
1246 if (ncontrib == 0) {
1247 outdata[ind] = fileptrs->sky;
1248 outdatac[ind] = 0;
1249 if (outvar != NULL)
1250 outdatav[ind] = 0.0;
1251 continue;
1252 }
1253 cc[ix] = 0;
1254
1255 /* Put some stuff in buffers for the averaging routine */
1256
1257 for (i = 0; i < ncontrib; i++) {
1258 dbuf[i] = workbuf[bufloc+i];
1259 wbuf[i] = workbufw[bufloc+i];
1260 cbuf[i] = workbufc[bufloc+i];
1261 if (outvar != NULL)
1262 vbuf[i] = workbufv[bufloc+i];
1263 ibuf[i] = id[bufloc+i];
1264 lbuf[i] = iloc[bufloc+i];
1265 }
1266
1267 /* Do the averages with the nominal clipping */
1268
1269 do_averages(ncontrib,dbuf,vbuf,wbuf,cbuf,ibuf,lclip,hclip,hsig,0.0,
1270 fileptrs,sumweights,&clipped_low,&clipped_high,&outval,
1271 &outvalc,&outvalv);
1272
1273 /* Store these away */
1274
1275 outdata[ind] = outval;
1276 confwork[ind] = outvalc;
1277 if (outvar != NULL)
1278 outdatav[ind] = outvalv;
1279 cc[ix] = (clipped_high >= 0);
1280 if (clipped_low >= 0)
1281 ((fileptrs+ibuf[clipped_low])->bpm)[lbuf[clipped_low]] = 2;
1282 if (clipped_high >= 0)
1283 ((fileptrs+ibuf[clipped_high])->bpm)[lbuf[clipped_high]] = 1;
1284 }
1285
1286 /* Look at the non-edge pixels and see if the local variation can
1287 show whether any clipping that has been done is justified. NB:
1288 we are looking at the _previous_ row here */
1289
1290 if (iy < 2)
1291 continue;
1292 for (ix = 1; ix < nxo-1; ix++) {
1293 if (! cclast[ix])
1294 continue;
1295 ind = (iy-1)*nxo + ix;
1296 bufloc = (iy-1)*nxo*nz + ix*nz;
1297 ncontrib = nbuf[ind];
1298 if (ncontrib == 0)
1299 continue;
1300 avlev = 0.0;
1301 for (i = -1; i <= 1; i++)
1302 for (j = -1; j <= 1; j++)
1303 avlev += outdata[ind + i*nxo + j];
1304 avlev /= 9.0;
1305 avvar = 0.0;
1306 for (i = -1; i <= 1; i++)
1307 for (j = -1; j <= 1; j++)
1308 avvar += (float)fabs(avlev - outdata[ind + i*nxo + j]);
1309 avvar *= scale_avvar;
1310 if (avlev <= hclip && avvar <= 2.0*(fileptrs->noise))
1311 continue;
1312
1313 /* Put some stuff in buffers for the averaging routine */
1314
1315 for (i = 0; i < ncontrib; i++) {
1316 dbuf[i] = workbuf[bufloc+i];
1317 if (outvar != NULL)
1318 vbuf[i] = workbufv[bufloc+i];
1319 wbuf[i] = workbufw[bufloc+i];
1320 cbuf[i] = workbufc[bufloc+i];
1321 ibuf[i] = id[bufloc+i];
1322 lbuf[i] = iloc[bufloc+i];
1323 (fileptrs+ibuf[i])->bpm[lbuf[i]] = 0;
1324 }
1325
1326 /* Redo the averages with the nominal clipping */
1327
1328 do_averages(ncontrib,dbuf,vbuf,wbuf,cbuf,ibuf,lclip,hclip,hsig,avvar,
1329 fileptrs,sumweights,&clipped_low,&clipped_high,&outval,
1330 &outvalc,&outvalv);
1331
1332 /* Store these away */
1333
1334 outdata[ind] = outval;
1335 if (outvar != NULL)
1336 outdatav[ind] = outvalv;
1337 confwork[ind] = outvalc;
1338 if (clipped_low >= 0)
1339 ((fileptrs+ibuf[clipped_low])->bpm)[lbuf[clipped_low]] = 2;
1340 if (clipped_high >= 0)
1341 ((fileptrs+ibuf[clipped_high])->bpm)[lbuf[clipped_high]] = 1;
1342 }
1343 }
1344
1345 /* Now renormalise the confidence map */
1346
1347 npts = nxo*nyo;
1348 bbpm = cpl_calloc(npts,sizeof(*bbpm));
1349 for (i = 0; i < npts; i++)
1350 bbpm[i] = (confwork[i] == 0.0);
1351 med = casu_med(confwork,bbpm,npts);
1352 freespace(bbpm);
1353 renorm = 100.0/med;
1354 for (i = 0; i < npts; i++) {
1355 confwork[i] *= renorm;
1356 outdatac[i] = max(0,min(1000,(int)(confwork[i]+0.5)));
1357 }
1358
1359 /* Free up some workspace */
1360
1361 freespace(workbuf);
1362 freespace(workbufc);
1363 freespace(workbufv);
1364 freespace(workbufw);
1365 freespace(confwork);
1366 freespace(id);
1367 freespace(iloc);
1368 freespace(nbuf);
1369 freespace(clipmap);
1370 freespace(dbuf);
1371 freespace(vbuf);
1372 freespace(wbuf);
1373 freespace(cbuf);
1374 freespace(ibuf);
1375 freespace(lbuf);
1376 return(CASU_OK);
1377}
1378
1379/*---------------------------------------------------------------------------*/
1409/*---------------------------------------------------------------------------*/
1410
1411static int stack_lin(int nim, dstrct *fileptrs, float sumweights,
1412 casu_fits *outfits, casu_fits *outconf, casu_fits *outvar,
1413 cpl_wcs *outwcs) {
1414 cpl_image *outim,*outimc,*outimv;
1415 double oxf,oyf;
1416 long npts,i,ix,iy,ind1,ind,ox,oy,m,jj,ii,nn,nx,ny,nxo,nyo;
1417 float *workbufw,*data,wt,dx,dy,w[4],value,cval,*vdata,*workdatav3;
1418 float med,*outdata,renorm,vval,*outdatav,*workdatav,*workdatav2;
1419 int *cdata,*outdatac,nin;
1420 dstrct *dd;
1421 unsigned char *bbpm;
1422 cpl_matrix *xyin,*xyout;
1423
1424 vdata = NULL; /* return value for cpl_image_get_data_float on error */
1425 workdatav = NULL; /* Remove falsepositive compiler warning from freespace macro*/
1426 workdatav2 = NULL; /* Remove falsepositive compiler warning from freespace macro*/
1427 workdatav3 = NULL; /* Remove falsepositive compiler warning from freespace macro*/
1428 outdatav = NULL;
1429 /* Output file info */
1430
1431 outim = casu_fits_get_image(outfits);
1432 outimc = casu_fits_get_image(outconf);
1433 if (outvar != NULL)
1434 outimv = casu_fits_get_image(outvar);
1435 nxo = (long)cpl_image_get_size_x(outim);
1436 nyo = (long)cpl_image_get_size_y(outim);
1437 npts = nxo*nyo;
1438
1439 /* Get workspace for collecting info over the whole output map */
1440
1441 workbufw = cpl_calloc(nxo*nyo,sizeof(*workbufw));
1442
1443 /* Get output buffers */
1444
1445 outdata = cpl_image_get_data_float(outim);
1446 outdatac = cpl_image_get_data_int(outimc);
1447 if (outvar != NULL) {
1448 workdatav = cpl_malloc(nxo*nyo*sizeof(float));
1449 workdatav2 = cpl_malloc(nxo*nyo*sizeof(float));
1450 workdatav3 = cpl_malloc(nxo*nyo*sizeof(float));
1451 outdatav = cpl_image_get_data_float(outimv);
1452 }
1453
1454 /* Initialise the output map. This is because the output data array
1455 probably still has the results of the nearest neighbour stack */
1456
1457 for (i = 0; i < npts; i++) {
1458 outdata[i] = 0.0;
1459 outdatac[i] = 0;
1460 if (outvar != NULL) {
1461 outdatav[i] = 0.0;
1462 workdatav3[i] = 0.0;
1463 }
1464 }
1465
1466 /* Ok, loop for each of the input files and collect all the information
1467 for each pixel */
1468
1469 for (i = 0; i < nim; i++) {
1470 if (outvar != NULL) {
1471 memset(workdatav,0,nxo*nyo*sizeof(float));
1472 memset(workdatav2,0,nxo*nyo*sizeof(float));
1473 }
1474 dd = fileptrs + i;
1475 data = cpl_image_get_data_float(casu_fits_get_image(dd->fname));
1476 cdata = cpl_image_get_data_int(casu_fits_get_image(dd->conf));
1477 if (outvar != NULL)
1478 vdata = cpl_image_get_data_float(casu_fits_get_image(dd->var));
1479 wt = dd->weight;
1480 nx = dd->nx;
1481 ny = dd->ny;
1482
1483 /* Get the position matrix */
1484
1485 nin = nx*ny;
1486 xyin = cpl_matrix_new((cpl_size)nin,2);
1487 nin = 0;
1488 for (iy = 1; iy <= ny; iy++) {
1489 for (ix = 1; ix <= nx; ix++) {
1490 cpl_matrix_set(xyin,nin,0,(double)ix);
1491 cpl_matrix_set(xyin,nin,1,(double)iy);
1492 nin++;
1493 }
1494 }
1495 outloc(dd->vwcs,xyin,outwcs,dd->trans,&xyout);
1496 cpl_matrix_delete(xyin);
1497
1498 /* Loop for each pixel in the data array and work out where it belongs
1499 in the output array */
1500
1501 nin = 0;
1502 for (iy = 0; iy < ny; iy++) {
1503 ind1 = iy*nx;
1504 for (ix = 0; ix < nx; ix++) {
1505 ind = ind1 + ix;
1506 oxf = cpl_matrix_get(xyout,nin,0);
1507 oyf = cpl_matrix_get(xyout,nin,1);
1508 nin++;
1509 if (oxf < 0.0 || oyf < 0.0)
1510 continue;
1511 ox = (long)oxf - 1;
1512 oy = (long)oyf - 1;
1513 if ((dd->bpm)[ind] != 0)
1514 continue;
1515 if (ox < 0 || ox >= nxo || oy < 0 || oy >= nyo ||
1516 cdata[ind] == 0)
1517 continue;
1518 dx = oxf - (float)(ox + 1);
1519 dy = oyf - (float)(oy + 1);
1520 w[0] = (1.0-dx)*(1.0-dy);
1521 w[1] = dx*(1.0-dy);
1522 w[2] = dy*(1.0-dx);
1523 w[3] = dx*dy;
1524 m = 0;
1525 value = data[ind]/(dd->expscale) + dd->skydiff;
1526 if (outvar != NULL)
1527 vval = vdata[ind]/((dd->expscale)*(dd->expscale));
1528 cval = (float)cdata[ind];
1529 for (jj = 0; jj <= 1; jj++) {
1530 for (ii = 0; ii <= 1; ii++) {
1531 nn = (oy+jj)*nxo + ox + ii;
1532 if (nn >= npts || ox+ii >= nxo || oy+jj >= nyo)
1533 continue;
1534 outdata[nn] += w[m]*value*wt*cval;
1535 workbufw[nn] += w[m]*wt*cval;
1536 if (outvar != NULL) {
1537 workdatav[nn] += vval*w[m]*cval;
1538 workdatav2[nn] += w[m]*cval;
1539 }
1540 m++;
1541 }
1542 }
1543 }
1544 }
1545 cpl_matrix_delete(xyout);
1546 if (outvar != NULL) {
1547 for (ii = 0; ii < npts; ii++) {
1548 if (workdatav2[ii] != 0.0) {
1549 outdatav[ii] += wt*wt*workdatav[ii]/workdatav2[ii];
1550 workdatav3[ii] += wt;
1551 }
1552 }
1553 }
1554 }
1555
1556 /* Now do the final averaging */
1557
1558 for (i = 0; i < npts; i++) {
1559 if (workbufw[i] != 0.0) {
1560 outdata[i] /= workbufw[i];
1561 if (outvar != NULL) {
1562 if (workdatav3[i] != 0.0)
1563 outdatav[i] /= powf(workdatav3[i],2.0);
1564 else
1565 outdatav[i] = 0.0;
1566 }
1567 workbufw[i] = workbufw[i]/sumweights;
1568 } else {
1569 outdata[i] = fileptrs->sky;
1570 if (outvar != NULL)
1571 outdatav[i] = 0.0;
1572 }
1573 }
1574
1575 /* Now renormalise the confidence map */
1576
1577 bbpm = cpl_calloc(npts,sizeof(*bbpm));
1578 for (i = 0; i < npts; i++)
1579 bbpm[i] = (workbufw[i] == 0.0);
1580 med = casu_med(workbufw,bbpm,npts);
1581 freespace(bbpm);
1582 renorm = 100.0/med;
1583 for (i = 0; i < npts; i++) {
1584 workbufw[i] *= renorm;
1585 outdatac[i] = max(0,min(1000,(int)(workbufw[i]+0.5)));
1586 }
1587 freespace(workbufw);
1588 if (outvar != NULL) {
1589 freespace(workdatav);
1590 freespace(workdatav2);
1591 freespace(workdatav3);
1592 }
1593 return(CASU_OK);
1594}
1595
1596/*---------------------------------------------------------------------------*/
1618/*---------------------------------------------------------------------------*/
1619
1620static void seeing_wt(int nim, dstrct *fileptrs) {
1621 float *seeing;
1622 int ierr,i;
1623 cpl_propertylist *plist;
1624 dstrct *dd;
1625 const char *fctid = "seeing_wt";
1626
1627 /* Get an array for the seeing estimates */
1628
1629 seeing = cpl_malloc(nim*sizeof(*seeing));
1630
1631 /* Read the seeing keyword. Get out of here if one is missing and
1632 just use the background weights that we already have */
1633
1634 ierr = 0;
1635 for (i = 0; i < nim; i++) {
1636 dd = fileptrs+i;
1637 plist = casu_fits_get_ehu(dd->fname);
1638 if (cpl_propertylist_has(plist,"ESO DRS SEEING")) {
1639 seeing[i] = cpl_propertylist_get_float(plist,"ESO DRS SEEING");
1640 if (seeing[i] <= 0.0) {
1641 ierr = 1;
1642 break;
1643 }
1644 } else {
1645 ierr = 1;
1646 break;
1647 }
1648 }
1649
1650 /* If there were no errors, then modify the background weights by the
1651 seeing ratio */
1652
1653 if (! ierr) {
1654 for (i = 1; i < nim; i++) {
1655 dd = fileptrs + i;
1656 dd->weight *= seeing[0]/seeing[i];
1657 }
1658 } else {
1659 cpl_msg_warning(fctid,"No seeing weighting will be done");
1660 }
1661 freespace(seeing);
1662}
1663
1664/*---------------------------------------------------------------------------*/
1690/*---------------------------------------------------------------------------*/
1691
1692static void output_files(int nim, dstrct *fileptrs, casu_fits **outfits,
1693 casu_fits **outconf, casu_fits **outvar,
1694 cpl_wcs **outwcs) {
1695 cpl_image *newim,*newconf,*newvar;
1696 double lowerleft[2],upperleft[2],lowerright[2],upperright[2];
1697 double xmin,ymin,sh[2],ra,dec;
1698 dstrct *dd;
1699 int i,ixo,iyo,isnull;
1700 cpl_wcs *refwcs;
1701 cpl_matrix *xyin,*xyout;
1702 cpl_propertylist *refp;
1703
1704 /* Initialise the first element of the results arrays since we're
1705 doing all this relative to the first image */
1706
1707 lowerleft[0] = 1.0;
1708 lowerleft[1] = 1.0;
1709 upperleft[0] = 1.0;
1710 upperleft[1] = (double)fileptrs->ny;
1711 lowerright[0] = (double)fileptrs->nx;
1712 lowerright[1] = 1.0;
1713 upperright[0] = (double)fileptrs->nx;
1714 upperright[1] = (double)fileptrs->ny;
1715
1716 /* Get a matrix to use for the coordinate transformation of the
1717 corner positions */
1718
1719 xyin = cpl_matrix_new(4,2);
1720
1721 /* Loop for all images relative to the first one */
1722
1723 refwcs = NULL;
1724 for (i = 1; i < nim; i++) {
1725 if (i == 1) {
1726 refp = cpl_propertylist_duplicate(casu_fits_get_ehu(fileptrs->fname));
1727 if (! strncmp(cpl_propertylist_get_string(refp,"CTYPE1"),"RA",2))
1728 cpl_propertylist_update_string(refp,"CTYPE1","RA---TAN");
1729 else
1730 cpl_propertylist_update_string(refp,"CTYPE2","RA---TAN");
1731 if (! strncmp(cpl_propertylist_get_string(refp,"CTYPE1"),"DEC",3))
1732 cpl_propertylist_update_string(refp,"CTYPE1","DEC--TAN");
1733 else
1734 cpl_propertylist_update_string(refp,"CTYPE2","DEC--TAN");
1735 refwcs = cpl_wcs_new_from_propertylist(refp);
1736 cpl_propertylist_delete(refp);
1737 }
1738 dd = fileptrs + i;
1739 cpl_matrix_set(xyin,0,0,(double)1.0);
1740 cpl_matrix_set(xyin,0,1,(double)1.0);
1741 cpl_matrix_set(xyin,1,0,(double)1.0);
1742 cpl_matrix_set(xyin,1,1,(double)(dd->ny));
1743 cpl_matrix_set(xyin,2,0,(double)(dd->nx));
1744 cpl_matrix_set(xyin,2,1,(double)1.0);
1745 cpl_matrix_set(xyin,3,0,(double)(dd->nx));
1746 cpl_matrix_set(xyin,3,1,(double)(dd->ny));
1747 outloc(dd->vwcs,xyin,refwcs,dd->trans,&xyout);
1748 lowerleft[0] = min(lowerleft[0],cpl_matrix_get(xyout,0,0));
1749 lowerleft[1] = min(lowerleft[1],cpl_matrix_get(xyout,0,1));
1750 upperleft[0] = min(upperleft[0],cpl_matrix_get(xyout,1,0));
1751 upperleft[1] = max(upperleft[1],cpl_matrix_get(xyout,1,1));
1752 lowerright[0] = max(lowerright[0],cpl_matrix_get(xyout,2,0));
1753 lowerright[1] = min(lowerright[1],cpl_matrix_get(xyout,2,1));
1754 upperright[0] = max(upperright[0],cpl_matrix_get(xyout,3,0));
1755 upperright[1] = max(upperright[1],cpl_matrix_get(xyout,3,1));
1756 cpl_matrix_delete(xyout);
1757 }
1758 freewcs(refwcs);
1759 cpl_matrix_delete(xyin);
1760
1761 /* Ok, what are the limits? */
1762
1763 ixo = casu_nint(max(lowerright[0]-lowerleft[0],upperright[0]-upperleft[0])) + 1;
1764 iyo = casu_nint(max(upperright[1]-lowerright[1],upperleft[1]-lowerleft[1])) + 1;
1765 xmin = min(lowerleft[0],upperleft[0]);
1766 ymin = min(lowerleft[1],lowerright[1]);
1767
1768 /* Create the image */
1769
1770 newim = cpl_image_new((cpl_size)ixo,(cpl_size)iyo,CPL_TYPE_FLOAT);
1771 *outfits = casu_fits_wrap(newim,fileptrs->fname,NULL,NULL);
1772 cpl_propertylist_update_int(casu_fits_get_ehu(*outfits),"NAXIS1",ixo);
1773 cpl_propertylist_update_int(casu_fits_get_ehu(*outfits),"NAXIS2",iyo);
1774
1775 /* Output file is always TAN projection */
1776
1777 refp = casu_fits_get_ehu(*outfits);
1778 cpl_propertylist_update_string(refp,"CTYPE1","RA---TAN");
1779 cpl_propertylist_update_string(refp,"CTYPE2","DEC--TAN");
1780 cpl_propertylist_erase_regexp(refp,"^PV[0-9]*_[0-9]*$",0);
1781
1782 /* Update the reference point for the WCS to take the jitter
1783 offsets into account */
1784
1785 sh[0] = xmin;
1786 sh[1] = ymin;
1787 (void)casu_crpixshift(refp,1.0,sh);
1788 *outwcs = cpl_wcs_new_from_propertylist(casu_fits_get_ehu(*outfits));
1789 ra = cpl_array_get(cpl_wcs_get_crval(*outwcs),0,&isnull);
1790 dec = cpl_array_get(cpl_wcs_get_crval(*outwcs),1,&isnull);
1791 casu_radectoxy(*outwcs,ra,dec,&xmin,&ymin);
1792 cpl_wcs_delete(*outwcs);
1793 casu_propertylist_update_float(refp,"CRPIX1",(float)xmin);
1794 casu_propertylist_update_float(refp,"CRPIX2",(float)ymin);
1795 *outwcs = cpl_wcs_new_from_propertylist(casu_fits_get_ehu(*outfits));
1796
1797 /* Now create the confidence map image */
1798
1799 newconf = cpl_image_new((cpl_size)ixo,(cpl_size)iyo,CPL_TYPE_INT);
1800 *outconf = casu_fits_wrap(newconf,*outfits,NULL,NULL);
1801
1802 /* Create the output variance map if required */
1803
1804 if (fileptrs->var != NULL) {
1805 newvar = cpl_image_new((cpl_size)ixo,(cpl_size)iyo,CPL_TYPE_FLOAT);
1806 *outvar = casu_fits_wrap(newvar,*outfits,NULL,NULL);
1807 } else {
1808 *outvar = NULL;
1809 }
1810}
1811
1812/*---------------------------------------------------------------------------*/
1835/*---------------------------------------------------------------------------*/
1836
1837static void diffxy(cpl_wcs *wref, cpl_wcs *wprog, float *dx, float *dy) {
1838 cpl_matrix *xyin,*xyout;
1839 double xref,yref;
1840 const cpl_array *arr;
1841 int isnull;
1842
1843 /* Find the middle of the reference frame */
1844
1845 arr = cpl_wcs_get_image_dims(wref);
1846 xref = 0.5*cpl_array_get(arr,0,&isnull);
1847 yref = 0.5*cpl_array_get(arr,1,&isnull);
1848 xyin = cpl_matrix_new(1,2);
1849 cpl_matrix_set(xyin,0,0,xref);
1850 cpl_matrix_set(xyin,0,1,yref);
1851
1852 /* Work out the output position */
1853
1854 outloc(wref,xyin,wprog,NULL,&xyout);
1855
1856 /* Now the offsets */
1857
1858 *dx = (float)(xref - cpl_matrix_get(xyout,0,0));
1859 *dy = (float)(yref - cpl_matrix_get(xyout,0,1));
1860 cpl_matrix_delete(xyout);
1861 cpl_matrix_delete(xyin);
1862}
1863
1864/*---------------------------------------------------------------------------*/
1889/*---------------------------------------------------------------------------*/
1890
1891static void outloc(cpl_wcs *win, cpl_matrix *in, cpl_wcs *wout,
1892 cpl_array *trans, cpl_matrix **out) {
1893 double *outdata,*tdata,xt,yt,xout,yout;
1894 cpl_size n;
1895 int i;
1896 cpl_array *status;
1897 cpl_matrix *radec;
1898
1899 /* First do the conversion of (x1,y1) => (ra,dec) */
1900
1901 cpl_wcs_convert(win,in,&radec,&status,CPL_WCS_PHYS2WORLD);
1902 cpl_array_delete(status);
1903
1904 /* Now do the conversion of (ra,dec) => (x2,y2) */
1905
1906 cpl_wcs_convert(wout,radec,out,&status,CPL_WCS_WORLD2PHYS);
1907 cpl_array_delete(status);
1908 cpl_matrix_delete(radec);
1909
1910 /* If there is a plate transformation, then we need to apply that
1911 here */
1912
1913 if (trans != NULL) {
1914 tdata = cpl_array_get_data_double(trans);
1915 n = cpl_matrix_get_nrow(*out);
1916 outdata = cpl_matrix_get_data(*out);
1917 for (i = 0; i < n; i++) {
1918 xt = outdata[2*i];
1919 yt = outdata[2*i+1];
1920 xout = xt*tdata[0] + yt*tdata[1] + tdata[2];
1921 yout = xt*tdata[3] + yt*tdata[4] + tdata[5];
1922 outdata[2*i] = xout;
1923 outdata[2*i+1] = yout;
1924 }
1925 }
1926}
1927
1928/*---------------------------------------------------------------------------*/
1953/*---------------------------------------------------------------------------*/
1954
1955static void inloc(cpl_wcs *wout, cpl_matrix *out, cpl_wcs *win,
1956 cpl_array *trans, cpl_matrix **in) {
1957 double *outdata,*tdata,xt,yt,xout,yout;
1958 cpl_size n;
1959 int i;
1960 cpl_matrix *radec;
1961 cpl_array *status;
1962
1963 /* If there is a plate transformation, then we need to apply that
1964 here */
1965
1966 if (trans != NULL) {
1967 tdata = cpl_array_get_data_double(trans);
1968 n = cpl_matrix_get_nrow(out);
1969 outdata = cpl_matrix_get_data(out);
1970 for (i = 0; i < n; i++) {
1971 xt = outdata[2*i];
1972 yt = outdata[2*i+1];
1973 xout = xt*tdata[0] + yt*tdata[1] + tdata[2];
1974 yout = xt*tdata[3] + yt*tdata[4] + tdata[5];
1975 outdata[2*i] = xout;
1976 outdata[2*i+1] = yout;
1977 }
1978 }
1979
1980 /* Now do the conversion of (x2,y2) => (ra,dec) */
1981
1982 cpl_wcs_convert(wout,out,&radec,&status,CPL_WCS_PHYS2WORLD);
1983 cpl_array_delete(status);
1984
1985 /* Now do the conversion of (ra,dec) => (x2,y2) */
1986
1987 cpl_wcs_convert(win,radec,in,&status,CPL_WCS_WORLD2PHYS);
1988 cpl_array_delete(status);
1989 cpl_matrix_delete(radec);
1990
1991}
1992
1993/*---------------------------------------------------------------------------*/
2015/*---------------------------------------------------------------------------*/
2016
2017static void backgrnd_ov(int nim, dstrct *fileptrs, int unmap) {
2018 int i,dx,dy,jj,nx,ii,*cdata1,*cdata,*cdata2;
2019 long npts1,npts,start[2],end[2],ind1,ind,n;
2020 dstrct *dd1,*dd2;
2021 float *data,sky1,skynoise1,sky2,skynoise2,frac,*data1,*data2;
2022
2023 /* Initialise a few things */
2024
2025 dd1 = fileptrs;
2026 npts1 = dd1->nx*dd1->ny;
2027
2028 /* Start by working out the full background of the first frame */
2029
2030 data1 = cpl_image_get_data_float(casu_fits_get_image(dd1->fname));
2031 cdata1 = cpl_image_get_data_int(casu_fits_get_image(dd1->conf));
2032 skyest(data1,cdata1,npts1,3.0,&sky1,&skynoise1);
2033 dd1->sky = sky1;
2034 dd1->noise = skynoise1;
2035 dd1->sky /= dd1->expscale;
2036 dd1->skydiff = 0.0;
2037 dd1->weight = 1.0;
2038
2039 /* Now loop for all the others */
2040
2041 for (i = 1; i < nim; i++) {
2042 dd2 = fileptrs + i;
2043 data2 = cpl_image_get_data_float(casu_fits_get_image(dd2->fname));
2044 cdata2 = cpl_image_get_data_int(casu_fits_get_image(dd2->conf));
2045
2046 /* Offset differences */
2047
2048 dx = casu_nint(dd1->xoff - dd2->xoff);
2049 dy = casu_nint(dd1->yoff - dd2->yoff);
2050
2051 /* Ok, here's the tricky bit. Work out the range of x and y for
2052 each of the images. The shift moves an image to the left and
2053 and down. Hence positive values of dx or dy moves the second
2054 frame further to the left of down compared to the first frame.
2055 Here are the coordinate ranges for the first image */
2056
2057 /* start[0] = max(0,-dx); */
2058 /* start[1] = max(0,-dy); */
2059 /* end[0] = min((dd1->nx-1),((dd1->nx)-dx-1)); */
2060 /* end[1] = min((dd1->ny-1),((dd1->ny)-dy-1)); */
2061 nx = dd1->nx;
2062 if (abs(dx) > dd1->nx || abs(dy) > dd1->ny) {
2063 frac = 0.0;
2064 } else {
2065 if (dx >= 0.0) {
2066 start[0] = 0;
2067 end[0] = max(0,dd1->nx-dx-1);
2068 } else {
2069 start[0] = abs(dx);
2070 end[0] = dd1->nx-1;
2071 }
2072 if (dy >= 0.0) {
2073 start[1] = 0;
2074 end[1] = max(0,dd1->ny-dy-1);
2075 } else {
2076 start[1] = abs(dy);
2077 end[1] = dd1->ny-1;
2078 }
2079
2080 /* How much does this cover the first image? */
2081
2082 npts = (end[0]-start[0]+1)*(end[1]-start[1]+1);
2083 frac = (float)npts/(float)npts1;
2084 }
2085
2086 /* If the coverage is more than 50% then calculate the background
2087 in the coverage region for the first image */
2088
2089 if (frac >= 0.5) {
2090 data = cpl_malloc(npts*sizeof(*data));
2091 cdata = cpl_malloc(npts*sizeof(*cdata));
2092 n = 0;
2093 for (jj = start[1]; jj <= end[1]; jj++) {
2094 ind1 = jj*nx;
2095 for (ii = start[0]; ii <= end[0]; ii++) {
2096 ind = ind1 + ii;
2097 data[n] = data1[ind];
2098 cdata[n++] = cdata1[ind];
2099 }
2100 }
2101 skyest(data,cdata,n,3.0,&sky1,&skynoise1);
2102 freespace(data);
2103 freespace(cdata);
2104 } else {
2105 sky1 = dd1->sky;
2106 }
2107
2108 /* And here are the coordinate ranges for the second image if
2109 the coverage is more than 50%. Get the subset you need. */
2110
2111 if (frac > 0.5) {
2112 start[0] = max(0,dx);
2113 start[1] = max(0,dy);
2114 end[0] = min((dd2->nx-1),((dd2->nx)+dx-1));
2115 end[1] = min((dd2->ny-1),((dd2->ny)+dy-1));
2116 nx = dd2->nx;
2117 npts = (end[0]-start[0]+1)*(end[1]-start[1]+1);
2118 data = cpl_malloc(npts*sizeof(*data));
2119 cdata = cpl_malloc(npts*sizeof(*cdata));
2120 n = 0;
2121 for (jj = start[1]; jj <= end[1]; jj++) {
2122 ind1 = jj*nx;
2123 for (ii = start[0]; ii <= end[0]; ii++) {
2124 ind = ind1 + ii;
2125 data[n] = data2[ind];
2126 cdata[n++] = cdata2[ind];
2127 }
2128 }
2129 skyest(data,cdata,n,3.0,&sky2,&skynoise2);
2130 freespace(data);
2131 freespace(cdata);
2132 } else {
2133 npts = (dd2->nx)*(dd2->ny);
2134 skyest(data2,cdata2,npts,3.0,&sky2,&skynoise2);
2135 }
2136
2137 /* Store info away */
2138
2139 dd2->sky = sky2/dd2->expscale;
2140 dd2->skydiff = sky1 - sky2;
2141 dd2->noise = skynoise2/dd2->expscale;
2142 if (dd2->noise != 0.0 && dd2->sky != 0.0)
2143 dd2->weight = (float)(pow(dd1->noise,2.0)/pow(dd2->noise,2.0));
2144 else
2145 dd2->weight = 0.0;
2146 if (unmap) {
2147 casu_fits_unload_im(dd2->fname);
2148 if (! dd2->confcopy)
2149 casu_fits_unload_im(dd2->conf);
2150 }
2151 }
2152 if (unmap) {
2153 casu_fits_unload_im(dd1->fname);
2154 casu_fits_unload_im(dd1->conf);
2155 }
2156}
2157
2158/*---------------------------------------------------------------------------*/
2206/*---------------------------------------------------------------------------*/
2207
2208static void do_averages(int ncontrib, float *data, float *vdata, float *wbuf,
2209 float *conf, unsigned char *id, float lclip, float hclip,
2210 float hsig, float extra, dstrct *fileptrs, float sumweights,
2211 short int *lowcl, short int *highcl, float *outval,
2212 float *outvalc, float *outvalv) {
2213 float minval,maxval,sum,wconfsum,reflev,clipval,clipscale,scl;
2214 float csum,vsum;
2215 int imax,imin,i;
2216
2217 /* Do the summation and keep track of the minimum and
2218 maximum values */
2219
2220 minval = 1.0e10;
2221 maxval = -1.0e10;
2222 sum = 0.0;
2223 vsum = 0.0;
2224 wconfsum = 0.0;
2225 csum = 0.0;
2226 imin = -1;
2227 imax = -1;
2228 for (i = 0; i < ncontrib; i++) {
2229 sum += data[i]*wbuf[i]*conf[i];
2230 wconfsum += wbuf[i]*conf[i];
2231 if (vdata != NULL)
2232 vsum += vdata[i]*wbuf[i]*wbuf[i]*conf[i]*conf[i];
2233 csum += conf[i];
2234 if (data[i] > maxval) {
2235 maxval = data[i];
2236 imax = i;
2237 }
2238 if (data[i] < minval) {
2239 minval = data[i];
2240 imin = i;
2241 }
2242 }
2243
2244 /* First crack at output value */
2245
2246 if (wconfsum > 0.0) {
2247 *outval = sum/wconfsum;
2248 } else {
2249 *outval = fileptrs->sky;
2250 }
2251
2252 /* Look at high clipping for of cosmic rays */
2253
2254 *highcl = -1;
2255 if (maxval > hclip && wconfsum > 150.0 && csum > 150.0 &&
2256 imin != -1 && imax != -1) {
2257 scl = wbuf[imax]*conf[imax];
2258 reflev = (sum - data[imax]*scl)/(wconfsum - scl);
2259 clipscale = max(1.0,reflev)/max(1.0,(fileptrs+id[imax])->sky);
2260 clipval = reflev + clipscale*hsig*(fileptrs+id[imax])->noise;
2261 clipval += hsig*extra;
2262
2263 /* Clip the maximum point */
2264
2265 if (maxval > clipval) {
2266 sum -= data[imax]*scl;
2267 wconfsum -= scl;
2268 if (vdata != NULL)
2269 vsum -= vdata[imax]*scl*scl;
2270 *outval = reflev;
2271 *highcl = imax;
2272 }
2273 }
2274
2275 /* Clip low points (presumably serious cosmetic issues on detectors) */
2276
2277 *lowcl = -1;
2278 if (minval < lclip && wconfsum > 150.0 && csum > 150.0 &&
2279 imin != -1 && imax != -1) {
2280 scl = wbuf[imin]*conf[imin];
2281 reflev = (sum - data[imin]*scl)/(wconfsum - scl);
2282 clipscale = max(1.0,reflev)/max(1.0,(fileptrs+id[imin])->sky);
2283 clipval = reflev + clipscale*hsig*(fileptrs+id[imin])->noise;
2284 if (minval < clipval) {
2285 sum -= data[imin]*scl;
2286 wconfsum -= scl;
2287 if (vdata != NULL)
2288 vsum -= vdata[imin]*scl*scl;
2289 *outval = reflev;
2290 *lowcl = imin;
2291 }
2292 }
2293
2294 /* Output confidence (not normalised) */
2295
2296 *outvalc = wconfsum/sumweights;
2297 if (vdata != NULL)
2298 *outvalv = vsum/powf(wconfsum,2.0);
2299 else
2300 *outvalv = 0.0;
2301}
2302
2303/*---------------------------------------------------------------------------*/
2329/*---------------------------------------------------------------------------*/
2330
2331static void skyest(float *data, int *cdata, long npts, float thresh,
2332 float *skymed, float *skynoise) {
2333 unsigned char *bpm;
2334 long i;
2335
2336 /* Get a dummy bad pixel mask */
2337
2338 bpm = cpl_calloc(npts,sizeof(*bpm));
2339 for (i = 0; i < npts; i++)
2340 bpm[i] = (cdata[i] == 0);
2341
2342 /* Get the stats */
2343
2344 casu_qmedsig(data,bpm,npts,thresh,2,DATAMIN,DATAMAX,skymed,skynoise);
2345
2346 /* Clean up */
2347
2348 freespace(bpm);
2349}
2350
2351/*---------------------------------------------------------------------------*/
2374/*---------------------------------------------------------------------------*/
2375
2376static cpl_table *mknewtab(casu_tfits *cat, dstrct *ddcur, dstrct *ddref) {
2377 cpl_table *newtab;
2378 int nr,i;
2379 float *x,*y;
2380 cpl_matrix *xyin,*xyout;
2381
2382 /* Make a copy of the input table */
2383
2384 newtab = cpl_table_duplicate(casu_tfits_get_table(cat));
2385 nr = (int)cpl_table_get_nrow(newtab);
2386
2387 /* Get the columns of interest */
2388
2389 x = cpl_table_get_data_float(newtab,"X_coordinate");
2390 y = cpl_table_get_data_float(newtab,"Y_coordinate");
2391 xyin = cpl_matrix_new((cpl_size)nr,2);
2392 for (i = 0; i < nr; i++) {
2393 cpl_matrix_set(xyin,i,0,x[i]);
2394 cpl_matrix_set(xyin,i,1,y[i]);
2395 }
2396
2397 /* Do the conversion */
2398
2399 outloc(ddcur->vwcs,xyin,ddref->vwcs,NULL,&xyout);
2400 cpl_matrix_delete(xyin);
2401
2402 /* Loop for each object */
2403
2404 for (i = 0; i < nr; i++) {
2405 x[i] = (float)cpl_matrix_get(xyout,i,0);
2406 y[i] = (float)cpl_matrix_get(xyout,i,1);
2407 }
2408 cpl_matrix_delete(xyout);
2409 return(newtab);
2410}
2411
2412/*---------------------------------------------------------------------------*/
2427/*---------------------------------------------------------------------------*/
2428
2429static cpl_array *transinit(void) {
2430 double *td;
2431 cpl_array *t;
2432
2433 t = cpl_array_new(6,CPL_TYPE_DOUBLE);
2434 td = cpl_array_get_data_double(t);
2435 td[0] = 1.0;
2436 td[1] = 0.0;
2437 td[2] = 0.0;
2438 td[3] = 0.0;
2439 td[4] = 1.0;
2440 td[5] = 0.0;
2441 return(t);
2442}
2443
2444/*---------------------------------------------------------------------------*/
2461/*---------------------------------------------------------------------------*/
2462
2463static void fileptrs_close(int nim, dstrct *fileptrs) {
2464 int i;
2465
2466 if (fileptrs == NULL)
2467 return;
2468 for (i = 0; i < nim; i++) {
2469 freewcs((fileptrs+i)->vwcs);
2470 freespace((fileptrs+i)->bpm);
2471 freearray((fileptrs+i)->trans);
2472 freearray((fileptrs+i)->invtrans);
2473 }
2474 cpl_free(fileptrs);
2475}
2476
2477static void tidy(int nim, dstrct *fileptrs, cpl_wcs *outwcs) {
2478 fileptrs_close(nim,fileptrs);
2479 if (outwcs != NULL)
2480 cpl_wcs_delete(outwcs);
2481}
2482
2483/*
2484
2485$Log: casu_imstack.c,v $
2486Revision 1.10 2015/11/25 10:27:19 jim
2487Modified to allow for unmapping
2488
2489Revision 1.9 2015/10/22 19:49:02 jim
2490Can't unload images just in they don't exist on disc to start with (because
2491then you can't reload them)
2492
2493Revision 1.8 2015/10/22 10:59:45 jim
2494Modified for speed and memory usage
2495
2496Revision 1.7 2015/09/14 18:46:26 jim
2497replaced a call to free with cpl_free
2498
2499Revision 1.6 2015/08/12 11:20:47 jim
2500removed some unnecessary variables
2501
2502Revision 1.5 2015/08/07 13:06:54 jim
2503Fixed copyright to ESO
2504
2505Revision 1.4 2015/08/06 05:34:02 jim
2506Fixes to get rid of compiler moans
2507
2508Revision 1.3 2015/08/03 12:44:05 jim
2509Fixed bug causing memory leak
2510
2511Revision 1.2 2015/06/25 11:53:27 jim
2512Fixed bug in overlap routine where the wrong section gets analysed
2513
2514Revision 1.1.1.1 2015/06/12 10:44:32 jim
2515Initial import
2516
2517Revision 1.24 2015/06/09 18:34:25 jim
2518Fixed bug in return value of variance for n=1 input
2519
2520Revision 1.23 2015/05/20 11:59:02 jim
2521Removed a redundant memory free
2522
2523Revision 1.22 2015/05/14 11:46:57 jim
2524Fixed bug in slow routine affecting the coordinate conversion
2525
2526Revision 1.21 2015/05/13 11:44:52 jim
2527Removed homemade WCS routines and replaced them with cpl standard ones
2528
2529Revision 1.20 2015/03/12 09:16:51 jim
2530Modified to remove some compiler moans
2531
2532Revision 1.19 2015/01/29 11:51:56 jim
2533modified comments
2534
2535Revision 1.18 2015/01/09 12:13:15 jim
2536*** empty log message ***
2537
2538Revision 1.17 2014/12/14 21:14:53 jim
2539Fixed a little bug
2540
2541Revision 1.16 2014/12/11 17:14:31 jim
2542Fixed issue when only 1 file is being stacked
2543
2544Revision 1.15 2014/12/11 12:23:33 jim
2545new version
2546
2547Revision 1.14 2014/04/26 18:41:22 jim
2548fixed lots of little bugs
2549
2550Revision 1.13 2014/04/24 10:26:51 jim
2551Fixed potential zero divide problem
2552
2553Revision 1.12 2014/04/23 05:40:46 jim
2554superficial modifications
2555
2556Revision 1.11 2014/04/09 11:08:21 jim
2557Get rid of a couple of compiler moans
2558
2559Revision 1.10 2014/04/09 09:09:51 jim
2560Detabbed
2561
2562Revision 1.9 2014/04/04 08:42:40 jim
2563Fixed scaling of output confidence
2564
2565Revision 1.8 2014/04/03 11:59:29 jim
2566Removed diagnostics
2567
2568Revision 1.7 2014/04/03 11:55:50 jim
2569Fixed scaling problem in rejection
2570
2571Revision 1.6 2014/03/26 15:48:07 jim
2572Modified to remove globals. Also augmented to put in a slower algorithm that
2573takes much less memory. First crack at doing error propagation.
2574
2575Revision 1.5 2013/11/21 09:38:14 jim
2576detabbed
2577
2578Revision 1.4 2013/11/21 08:51:35 jim
2579Removed some calls to casu_fits_unload_im
2580
2581Revision 1.3 2013/11/08 06:53:41 jim
2582Modified to try and minimise the amount of memory required
2583
2584Revision 1.2 2013-10-24 09:21:18 jim
2585modified the way the WCS is done
2586
2587Revision 1.1.1.1 2013-08-27 12:07:48 jim
2588Imported
2589
2590
2591
2592*/
casu_fits * casu_fits_wrap(cpl_image *im, casu_fits *model, cpl_propertylist *phu, cpl_propertylist *ehu)
Definition: casu_fits.c:883
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
casu_fits * casu_fits_duplicate(casu_fits *in)
Definition: casu_fits.c:225
cpl_propertylist * casu_fits_get_phu(casu_fits *p)
Definition: casu_fits.c:531
char * casu_fits_get_filename(casu_fits *p)
Definition: casu_fits.c:646
void casu_fits_unload_im(casu_fits *in)
Definition: casu_fits.c:200
cpl_propertylist * casu_fits_get_ehu(casu_fits *p)
Definition: casu_fits.c:576
int casu_fits_get_nexten(casu_fits *p)
Definition: casu_fits.c:497
int casu_platexy(cpl_table *matchedxy, int nconst, cpl_array **coefs, int *status)
Work out an xy transformation for two coodinate lists.
int casu_imstack(casu_fits **inf, casu_fits **inconf, casu_fits **invar, casu_tfits **cats, int nimages, int nconfs, float lthr, float hthr, int method, int seeing, int fast, int unmap, const char *expkey, casu_fits **out, casu_fits **outc, casu_fits **outv, int *status)
Stack images into a mean image using WCS info.
Definition: casu_imstack.c:187
int casu_matchxy(cpl_table *progtab, cpl_table *template, float srad, float *xoffset, float *yoffset, int *nm, cpl_table **outtab, int *status)
Match two lists of x,y coordinates from two tables to find the cartesian offset between them.
Definition: casu_match.c:105
float casu_med(float *data, unsigned char *bpm, long npts)
Definition: casu_stats.c:89
void casu_qmedsig(float *data, unsigned char *bpm, long npts, float thresh, int niter, float lowv, float highv, float *median, float *sigma)
Definition: casu_stats.c:258
cpl_table * casu_tfits_get_table(casu_tfits *p)
Definition: casu_tfits.c:364
void casu_prov(cpl_propertylist *p, casu_fits **inlist, int n, int isextn)
Write provenance keywords.
Definition: casu_utils.c:287
void casu_propertylist_update_float(cpl_propertylist *plist, const char *name, float val)
Definition: casu_utils.c:1163
void casu_radectoxy(cpl_wcs *wcs, double ra, double dec, double *x, double *y)
int casu_crpixshift(cpl_propertylist *p, double scalefac, double sh[])
int casu_coverage(cpl_propertylist *plist, int fudge, double *ra1, double *ra2, double *dec1, double *dec2, int *status)