VIRCAM Pipeline 2.3.12
casu_nebuliser.c
1/* $Id: casu_nebuliser.c,v 1.3 2015/11/25 10:26:31 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:26:31 $
24 * $Revision: 1.3 $
25 * $Name: $
26 */
27
28/* Includes */
29
30#ifdef HAVE_CONFIG_H
31#include <config.h>
32#endif
33
34#include <cpl.h>
35#include "casu_mods.h"
36#include "catalogue/casu_utils.h"
37#include "catalogue/casu_fits.h"
38#include "casu_stats.h"
39
40#include <stdio.h>
41#include <stdlib.h>
42#include <string.h>
43#include <math.h>
44
45#define MEDIANCALC 1
46#define MEANCALC 2
47#define DATAMIN -1000.0
48#define DATAMAX 65535.0
49
50typedef struct {
51 float sum;
52 float sumw;
53 int naver;
54 float nextw;
55 float lastw;
56 float nextval;
57 float lastval;
58 short int nextc;
59 short int lastc;
60} nextlast;
61
62static void twodfilt(float *data, unsigned char *bpm, int nx, int ny,
63 int medfilt, int linfilt, int niter, int axis,
64 int twod, int takeout_sky, int inorm, int wantback,
65 float signeg, float sigpos, float **backmap);
66static void bfilt_1d(float *data, unsigned char *bpm, int nx, int ny,
67 int filt, int stat, int axis);
68static void bfilt_2d(float *data, unsigned char *bpmcopy, int nx, int ny,
69 int filt, int stat);
70static void docols_2(float *data, unsigned char *bpm, float *dbuf,
71 unsigned char *bbuf, int nx, int ny, int filter,
72 int stat);
73static void dorows_2(float *data, unsigned char *bpm, float *dbuf,
74 unsigned char *bbuf, int nx, int ny, int filter,
75 int stat);
76static void dostat(float *data, unsigned char *bpm, unsigned char *goodval,
77 int npts, int nfilt, int whichstat);
78static void wraparound(float *data, unsigned char *bpm, int npts, int nfilt,
79 int whichstat, float **ybuf, unsigned char **ybbuf,
80 int *nbuf);
81static void medavg(float *array, unsigned char *bpm, int *ipoint, int npix,
82 int whichstat, int newl, nextlast *nl, float *outval,
83 unsigned char *outbp);
84static void quickie(float *array, unsigned char *iarray, int *iarray2,
85 int lll, int narray);
86static void sortm(float *a1, unsigned char *a2, int *a3, int n);
87static void plugholes(float *data, unsigned char *bpm, int nx);
88
89
92/*---------------------------------------------------------------------------*/
151/*---------------------------------------------------------------------------*/
152
153extern int casu_nebuliser(casu_fits *infile, casu_fits *inconf, int medfilt,
154 int linfilt, int niter, int axis, int twod,
155 int takeout_sky, int norm, int wantback,
156 float signeg, float sigpos, casu_fits **backmap,
157 int *status) {
158 cpl_image *im;
159 int nx,ny,i,*cdata;
160 float *data,*backdata,*bdata;
161 int npts;
162 unsigned char *bpm;
163 const char *fctid = "casu_nebuliser";
164
165 /* Inherited status */
166
167 *backmap = NULL;
168 if (*status != CASU_OK)
169 return(*status);
170
171 /* Get the data array of the input file */
172
173 im = casu_fits_get_image(infile);
174 nx = cpl_image_get_size_x(im);
175 ny = cpl_image_get_size_y(im);
176 data = cpl_image_get_data_float(im);
177
178 /* Do we have a confidence map? If so then use it to create a bpm.
179 If not, then create a bpm with all zeros */
180
181 npts = nx*ny;
182 bpm = cpl_calloc(npts,sizeof(unsigned char));
183 if (inconf != NULL) {
184 im = casu_fits_get_image(inconf);
185 if (cpl_image_get_size_x(im) != nx || cpl_image_get_size_y(im) != ny) {
186 cpl_msg_error(fctid,"Image and conf map dimensions don't match");
187 freespace(bpm);
188 FATAL_ERROR
189 }
190 cdata = cpl_image_get_data(im);
191 for (i = 0; i < npts; i++)
192 bpm[i] = (cdata[i] == 0);
193 }
194
195 /* Ok do the filtering now */
196
197 twodfilt(data,bpm,nx,ny,medfilt,linfilt,niter,axis,twod,takeout_sky,norm,
198 wantback,signeg,sigpos,&backdata);
199
200 /* Write a header item to the data header to show that the nebulisation
201 has been successful */
202
203 cpl_propertylist_append_bool(casu_fits_get_ehu(infile),
204 "ESO DRS NEBULISED",1);
205 cpl_propertylist_set_comment(casu_fits_get_ehu(infile),
206 "ESO DRS NEBULISED",
207 "Nebuliser has been used on this image");
208
209 /* Now if we have a background image, then create the fits structure
210 for it as a copy of the input image */
211
212 if (wantback) {
213 *backmap = casu_fits_duplicate(infile);
214 im = casu_fits_get_image(*backmap);
215 bdata = cpl_image_get_data_float(im);
216 memmove(bdata,backdata,npts*sizeof(float));
217 freespace(backdata);
218 }
219
220 /* Tidy and exit */
221
222 freespace(bpm);
223 GOOD_STATUS
224}
225
226/*---------------------------------------------------------------------------*/
272/*---------------------------------------------------------------------------*/
273
274static void twodfilt(float *data, unsigned char *bpm, int nx, int ny,
275 int medfilt, int linfilt, int niter, int axis,
276 int twod, int takeout_sky, int inorm, int wantback,
277 float signeg, float sigpos, float **backmap) {
278 int i,iter,nter,nmask;
279 long nn;
280 float *buffer,*orig,*orig_sm,*work,medsky,sigsky,rescale,lthr,hthr;
281 float diff;
282 unsigned char *bpmcopy;
283
284 /* Get some workspace. One holds a copy of the original data. The
285 others are for work */
286
287 nn = nx*ny;
288 buffer = cpl_malloc(3*nn*sizeof(*buffer));
289 orig = buffer;
290 orig_sm = orig + nn;
291 work = orig_sm + nn;
292 memmove((char *)orig,(char *)data,nn*sizeof(*data));
293 memmove((char *)orig_sm,(char *)data,nn*sizeof(*data));
294 memmove((char *)work,(char *)data,nn*sizeof(*data));
295
296 /* Copy the bad pixel mask, so that the pre-existing bad pixels are
297 now flagged with 1. */
298
299 bpmcopy = cpl_calloc(nn,sizeof(*bpmcopy));
300 for (i = 0; i < nn; i++)
301 bpmcopy[i] = (bpm[i] ? 1 : 0);
302
303 /* Do gentle smooth on the original data */
304
305 if (niter > 1 && medfilt > 10) {
306 bfilt_1d(orig_sm,bpmcopy,nx,ny,5,MEDIANCALC,axis);
307 bfilt_1d(orig_sm,bpmcopy,nx,ny,3,MEANCALC,axis);
308 }
309
310 /* Now do an iteration loop */
311
312 for (iter = 1; iter <= niter; iter++) {
313 if (iter > 1)
314 memmove((char *)data,(char *)orig,nn*sizeof(*data));
315
316 /* Filter the original input data, using the latest interation
317 on the pixel mask */
318
319 if (! twod)
320 bfilt_1d(data,bpmcopy,nx,ny,medfilt,MEDIANCALC,axis);
321 else
322 bfilt_2d(data,bpmcopy,nx,ny,medfilt,MEDIANCALC);
323 if (iter == niter)
324 break;
325
326 /* Look at the difference between the smoothed map and the (possibly
327 gently smoothed) original data */
328
329 for (i = 0; i < nn; i++)
330 work[i] = orig_sm[i] - data[i];
331
332 /* What is the median level and RMS of the residual map? We may need
333 to iterate on this */
334
335 casu_qmedsig(work,bpmcopy,nn,3.0,3,DATAMIN,DATAMAX,&medsky,&sigsky);
336 rescale = 2.0;
337 nter = 0;
338 while (sigsky < 2.5 && nter < 16) {
339 nter++;
340 for (i = 0; i < nn; i++)
341 work[i] *= rescale;
342 casu_qmedsig(work,bpmcopy,nn,3.0,3,DATAMIN,DATAMAX,&medsky,
343 &sigsky);
344 }
345 if (nter > 0) {
346 rescale = (float)pow(2.0,(double)nter);
347 for (i = 0; i < nn; i++)
348 work[i] /= rescale;
349 medsky /= rescale;
350 sigsky /= rescale;
351 }
352 lthr = -signeg*sigsky;
353 hthr = sigpos*sigsky;
354
355 /* Clip out discordant points */
356
357 nmask = 0;
358 for (i = 0; i < nn; i++) {
359 if (bpmcopy[i] == 1)
360 continue;
361 diff = work[i] - medsky;
362 if (diff > hthr || diff < lthr) {
363 bpmcopy[i] = 2;
364 nmask++;
365 } else {
366 bpmcopy[i] = 0;
367 }
368 }
369 }
370
371 /* Now do the linear filter */
372
373 if (! twod)
374 bfilt_1d(data,bpm,nx,ny,linfilt,MEANCALC,axis);
375 else
376 bfilt_2d(data,bpm,nx,ny,linfilt,MEANCALC);
377
378 /* Get the sky level if you want to keep it */
379
380 if (! takeout_sky)
381 casu_qmedsig(orig,bpmcopy,nn,3.0,3,DATAMIN,DATAMAX,&medsky,&sigsky);
382 else
383 medsky = 0.0;
384
385 /* Do we want a background map */
386
387 if (wantback) {
388 *backmap = cpl_malloc(nn*sizeof(**backmap));
389 for (i = 0; i < nn; i++)
390 (*backmap)[i] = data[i];
391 } else {
392 *backmap = NULL;
393 }
394
395 /* How do we want to normalise? */
396
397 if (inorm == 0) {
398 for (i = 0; i < nn; i++)
399 data[i] = orig[i] - data[i] + medsky;
400 } else {
401 for (i = 0; i < nn; i++)
402 data[i] = orig[i]/max(1.0,data[i]);
403 }
404
405 /* Tidy and exit */
406
407 freespace(buffer);
408 freespace(bpmcopy);
409}
410
411/*---------------------------------------------------------------------------*/
439/*---------------------------------------------------------------------------*/
440
441static void bfilt_1d(float *data, unsigned char *bpm, int nx, int ny,
442 int filt, int stat, int axis) {
443
444 float *dbuf;
445 unsigned char *bbuf;
446 int nbuf;
447
448 /* Get some workspace */
449
450 nbuf = max(nx,ny);
451 dbuf = cpl_malloc(nbuf*sizeof(*dbuf));
452 bbuf = cpl_malloc(nbuf*sizeof(*bbuf));
453
454 /* Order the reset correction so that the first smoothing is done
455 across the axis of the anomaly */
456
457 if (axis == 1) {
458 dorows_2(data,bpm,dbuf,bbuf,nx,ny,filt,stat);
459 docols_2(data,bpm,dbuf,bbuf,nx,ny,filt,stat);
460 } else {
461 docols_2(data,bpm,dbuf,bbuf,nx,ny,filt,stat);
462 dorows_2(data,bpm,dbuf,bbuf,nx,ny,filt,stat);
463 }
464
465 /* Ditch workspace */
466
467 freespace(dbuf);
468 freespace(bbuf);
469}
470
471/*---------------------------------------------------------------------------*/
496/*---------------------------------------------------------------------------*/
497
498static void bfilt_2d(float *data, unsigned char *bpmcopy, int nx, int ny,
499 int filt, int stat) {
500 float *dbuf,*outmap,value,*om;
501 unsigned char *outbpm,*ob;
502 int nbuf,j,i,nf2,nalloc,jj,ii,ind,ind1,j1old,j2old,i1old,i2old;
503
504 /* Filter halfwidth */
505
506 nf2 = filt/2;
507
508 /* Get some workspace */
509
510 nalloc = (2*filt+1)*(2*filt+1);
511 dbuf = cpl_malloc(nalloc*sizeof(*dbuf));
512 outmap = cpl_malloc(nx*ny*sizeof(*outmap));
513 outbpm = cpl_malloc(nx*ny*sizeof(*outbpm));
514
515 /* Loop for each input pixel */
516
517 for (j = 0; j < ny; j++) {
518 for (i = 0; i < nx; i++) {
519 nbuf = 0;
520 for (jj = j - nf2; jj <= j + nf2; jj++) {
521 if (jj < 0 || jj >= ny)
522 continue;
523 ind1 = jj*nx;
524 for (ii = i - nf2; ii <= i + nf2; ii++) {
525 if (ii < 0 || ii >= nx)
526 continue;
527 ind = ind1 + ii;
528 if (bpmcopy[ind])
529 continue;
530 dbuf[nbuf++] = data[ind];
531 }
532 }
533
534 /* If we don't have enough, try and increase the window size.
535 This will only affect the edges */
536
537 if (nbuf < filt/4) {
538 j1old = j - nf2;
539 j2old = j + nf2;
540 i1old = i - nf2;
541 i2old = i + nf2;
542 for (jj = j - filt; jj <= j + filt; jj++) {
543 if (jj < 0 || jj >= ny || (jj >= j1old && jj <= j2old))
544 continue;
545 ind1 = jj*nx;
546 for (ii = i - filt; ii <= i + filt; ii++) {
547 if (ii < 0 || ii >= nx || (ii >= i1old && ii <= i2old))
548 continue;
549 ind = ind1 + ii;
550 if (bpmcopy[ind])
551 continue;
552 dbuf[nbuf++] = data[ind];
553 }
554
555 }
556 }
557
558 /* Right, assuming we have enough entries, then get a median */
559
560 ind = j*nx + i;
561 if (nbuf > filt/4) {
562 if (stat == MEDIANCALC)
563 value = casu_med(dbuf,NULL,(long)nbuf);
564 else
565 value = casu_mean(dbuf,NULL,(long)nbuf);
566 outmap[ind] = value;
567 outbpm[ind] = 0;
568 } else {
569 outmap[ind] = -1000.0;
570 outbpm[ind] = 1;
571 }
572 }
573 }
574
575 /* Right, fill in the holes and then transfer the filtered data to
576 the input array */
577
578 for (j = 0; j < ny; j++) {
579 om = outmap + j*nx;
580 ob = outbpm + j*nx;
581 plugholes(om,ob,nx);
582 for (i = 0; i < nx; i++)
583 data[j*nx+i] = om[i];
584 }
585
586 /* Tidy up and get out of here */
587
588 freespace(outmap);
589 freespace(outbpm);
590 freespace(dbuf);
591 return;
592}
593
594/*---------------------------------------------------------------------------*/
625/*---------------------------------------------------------------------------*/
626
627static void docols_2(float *data, unsigned char *bpm, float *dbuf,
628 unsigned char *bbuf, int nx, int ny, int filter,
629 int stat) {
630
631 int j,k,indx,nn;
632 unsigned char *goodval,*b;
633 float *t;
634
635 if (filter <= 0)
636 return;
637
638 goodval = cpl_malloc(ny*sizeof(*goodval));
639 t = cpl_malloc(ny*sizeof(*t));
640 b = cpl_malloc(ny*sizeof(*b));
641 for (k = 0; k < nx; k++) {
642 memset((char *)goodval,0,ny);
643 nn = 0;
644 for (j = 0; j < ny; j++) {
645 indx = j*nx + k;
646 if (bpm[indx] == 0) {
647 dbuf[nn] = data[indx];
648 bbuf[nn++] = 0;
649 }
650 }
651 dostat(dbuf,bbuf,goodval,nn,filter,stat);
652 nn = 0;
653 for (j = 0; j < ny; j++) {
654 indx = j*nx + k;
655 if (bpm[indx] == 0) {
656 t[j] = dbuf[nn++];
657 b[j] = 0;
658 } else {
659 t[j] = -999.0;
660 b[j] = 1;
661 }
662 }
663 plugholes(t,b,ny);
664 nn = 0;
665 for (j = 0; j < ny; j++) {
666 indx = j*nx + k;
667 data[indx] = t[j];
668 }
669 }
670 freespace(goodval);
671 freespace(t);
672 freespace(b);
673}
674
675/*---------------------------------------------------------------------------*/
706/*---------------------------------------------------------------------------*/
707
708static void dorows_2(float *data, unsigned char *bpm, float *dbuf,
709 unsigned char *bbuf, int nx, int ny, int filter,
710 int stat) {
711
712 int j,k,indx,nn;
713 unsigned char *goodval,*b;
714 float *t;
715
716 if (filter <= 0)
717 return;
718
719 goodval = cpl_malloc(nx*sizeof(*goodval));
720 t = cpl_malloc(nx*sizeof(*t));
721 b = cpl_malloc(nx*sizeof(*b));
722 for (k = 0; k < ny; k++) {
723 memset((char *)goodval,0,nx);
724 nn = 0;
725 for (j = 0; j < nx; j++) {
726 indx = k*nx + j;
727 if (bpm[indx])
728 continue;
729 dbuf[nn] = data[indx];
730 bbuf[nn++] = 0;
731 }
732 dostat(dbuf,bbuf,goodval,nn,filter,stat);
733 nn = 0;
734 for (j = 0; j < nx; j++) {
735 indx = k*nx + j;
736 if (bpm[indx] == 0) {
737 t[j] = dbuf[nn++];
738 b[j] = 0;
739 } else {
740 t[j] = -999.0;
741 b[j] = 1;
742 }
743 }
744 plugholes(t,b,nx);
745 for (j = 0; j < nx; j++) {
746 indx = k*nx + j;
747 data[indx] = t[j];
748 }
749 }
750 freespace(goodval);
751 freespace(t);
752 freespace(b);
753}
754
755/*---------------------------------------------------------------------------*/
783/*---------------------------------------------------------------------------*/
784
785static void dostat(float *data, unsigned char *bpm, unsigned char *goodval,
786 int npts, int nfilt, int whichstat) {
787 int nbuf,jl,jh,j,*ipoint,ifree,i;
788 unsigned char *ybbuf,*barray,bval;
789 nextlast nl;
790 float *ybuf,*darray,val;
791
792 /* Check to make sure there are some points in the array... */
793
794 if (npts < nfilt || npts < 10)
795 return;
796
797 /* Check to make sure the filter size is odd */
798
799 if ((nfilt/2)*2 == nfilt)
800 nfilt++;
801
802 /* Do the wrap around and load the data into an oversized array*/
803
804 wraparound(data,bpm,npts,nfilt,whichstat,&ybuf,&ybbuf,&nbuf);
805
806 /* Start doing the filtering...*/
807
808 darray = cpl_malloc(nfilt*sizeof(*darray));
809 barray = cpl_malloc(nfilt*sizeof(*barray));
810 ipoint = cpl_malloc(nfilt*sizeof(*ipoint));
811 memmove((char *)darray,(char *)ybuf,nfilt*sizeof(*ybuf));
812 memmove((char *)barray,(char *)ybbuf,nfilt*sizeof(*ybbuf));
813 for (j = 0; j < nfilt; j++)
814 ipoint[j] = j;
815 ifree = 0;
816 medavg(darray,barray,ipoint,nfilt,whichstat,-1,&nl,&val,&bval);
817 if (! bval)
818 data[0] = val;
819 goodval[0] = bval;
820 jl = nfilt;
821 jh = nfilt + npts - 2;
822 for (j = jl; j <= jh; j++) {
823 for (i = 0; i < nfilt; i++) {
824 if (ipoint[i] == 0) {
825 ifree = i;
826 ipoint[i] = nfilt - 1;
827 nl.lastval = darray[ifree];
828 nl.lastw = 0.0;
829 nl.lastc = 0;
830 if (barray[ifree] == 0) {
831 nl.lastw = 1.0;
832 nl.lastc = 1;
833 }
834 darray[ifree] = ybuf[j];
835 barray[ifree] = ybbuf[j];
836 nl.nextval = darray[ifree];
837 nl.nextw = 0.0;
838 nl.nextc = 0;
839 if (barray[ifree] == 0) {
840 nl.nextw = 1.0;
841 nl.nextc = 1;
842 }
843 } else
844 ipoint[i]--;
845 }
846 medavg(darray,barray,ipoint,nfilt,whichstat,ifree,&nl,&val,&bval);
847 if (! bval)
848 data[j-jl+1] = val;
849 goodval[j-jl+1] = bval;
850 }
851
852 /* Ditch workspace */
853
854 freespace(darray);
855 freespace(barray);
856 freespace(ipoint);
857 freespace(ybuf);
858 freespace(ybbuf);
859}
860
861/*---------------------------------------------------------------------------*/
891/*---------------------------------------------------------------------------*/
892
893static void wraparound(float *data, unsigned char *bpm, int npts, int nfilt,
894 int whichstat, float **ybuf, unsigned char **ybbuf,
895 int *nbuf) {
896
897 float *darray,xmns,xmnf;
898 int i1,ilow,i,*ipoint;
899 unsigned char *barray,bxmns,bxmnf;
900 nextlast nl;
901
902 /* Do some padding at the edges */
903
904 i1 = nfilt/2;
905 ilow = max(3,nfilt/4);
906 ilow = (ilow/2)*2 + 1;
907
908 /* Get some workspace */
909
910 darray = cpl_malloc(nfilt*sizeof(*darray));
911 barray = cpl_malloc(nfilt*sizeof(*barray));
912 ipoint = cpl_calloc(nfilt,sizeof(*ipoint));
913 *nbuf = npts + 2*i1;
914 *ybuf = cpl_malloc(*nbuf*sizeof(float));
915 *ybbuf = cpl_malloc(*nbuf*sizeof(unsigned char));
916
917 /* Do the wrap around.*/
918
919 memmove((char *)darray,(char *)data,ilow*sizeof(*data));
920 memmove((char *)barray,(char *)bpm,ilow*sizeof(*bpm));
921 medavg(darray,barray,ipoint,ilow,whichstat,-1,&nl,&xmns,&bxmns);
922 memmove((char *)darray,(char *)(data+npts-ilow),ilow*sizeof(*data));
923 memmove((char *)barray,(char *)(bpm+npts-ilow),ilow*sizeof(*bpm));
924 medavg(darray,barray,ipoint,ilow,whichstat,-1,&nl,&xmnf,&bxmnf);
925 for (i = 0; i < i1; i++) {
926 if (! bxmns) {
927 (*ybuf)[i] = 2.0*xmns - data[i1+ilow-i-1];
928 (*ybbuf)[i] = bpm[i1+ilow-i-1];
929 } else {
930 (*ybuf)[i] = data[i1+ilow-i-1];
931 (*ybbuf)[i] = 1;
932 }
933 if (! bxmnf) {
934 (*ybuf)[npts+i1+i] = 2.0*xmnf - data[npts-i-ilow-1];
935 (*ybbuf)[npts+i1+i] = bpm[npts-i-ilow-1];
936 } else {
937 (*ybuf)[npts+i1+i] = data[npts-i-ilow-1];
938 (*ybbuf)[npts+i1+i] = 1;
939 }
940 }
941
942 /* Now place the full line into the buffer */
943
944 memmove((char *)(*ybuf+i1),data,npts*sizeof(*data));
945 memmove((char *)(*ybbuf+i1),bpm,npts*sizeof(*bpm));
946
947 /* Free workspace */
948
949 freespace(darray);
950 freespace(barray);
951 freespace(ipoint);
952}
953
954/*---------------------------------------------------------------------------*/
986/*---------------------------------------------------------------------------*/
987
988static void medavg(float *array, unsigned char *bpm, int *ipoint, int npix,
989 int whichstat, int newl, nextlast *nl, float *outval,
990 unsigned char *outbp) {
991 float *buf = NULL;
992 int m,i;
993
994 /* If there is a new element and that new element is bad, then
995 the buffer is already sorted from last time (just one element
996 sorter */
997
998 m = 0;
999 if (whichstat == MEDIANCALC) {
1000 if (newl == -1)
1001 sortm(array,bpm,ipoint,npix);
1002 else
1003 quickie(array,bpm,ipoint,newl,npix);
1004
1005 /* Get some workspace */
1006
1007 buf = cpl_malloc(npix*sizeof(*buf));
1008
1009 /* Now put everything that's good in the buffer */
1010
1011 m = 0;
1012 for (i = 0; i < npix; i++) {
1013 if (bpm[i] == 0) {
1014 buf[m] = array[i];
1015 m++;
1016 }
1017 }
1018 } else if (whichstat == MEANCALC) {
1019 if (newl == -1) {
1020 nl->sum = 0.0;
1021 nl->sumw = 0.0;
1022 nl->naver = 0;
1023 for (i = 0; i < npix; i++) {
1024 if (bpm[i] == 0) {
1025 nl->sum += array[i];
1026 nl->sumw += 1.0;
1027 nl->naver += 1;
1028 }
1029 }
1030 m = nl->naver;
1031 } else {
1032 nl->sum += (nl->nextw*nl->nextval - nl->lastw*nl->lastval);
1033 nl->sumw += (nl->nextw - nl->lastw);
1034 nl->naver += (nl->nextc - nl->lastc);
1035 m = nl->naver;
1036 }
1037 }
1038
1039 /* If they were all bad, then send a null result back */
1040
1041 if (m == 0) {
1042 *outval = 0.0;
1043 *outbp = 1;
1044 if (whichstat == MEDIANCALC)
1045 freespace(buf);
1046
1047 /* Otherwise calculate the relevant stat */
1048
1049 } else {
1050 if (whichstat == MEDIANCALC) {
1051 *outval = buf[m/2];
1052 freespace(buf);
1053 } else if (whichstat == MEANCALC)
1054 *outval = (nl->sum)/(nl->sumw);
1055 *outbp = 0;
1056 }
1057}
1058
1059/*---------------------------------------------------------------------------*/
1085/*---------------------------------------------------------------------------*/
1086
1087static void quickie(float *array, unsigned char *iarray, int *iarray2,
1088 int lll, int narray) {
1089
1090 float test;
1091 int i,j,npt,it2;
1092 unsigned char it;
1093
1094 test = array[lll];
1095 it = iarray[lll];
1096 it2 = iarray2[lll];
1097 j = -1;
1098 for (i = 0; i < narray; i++) {
1099 if (i != lll && test <= array[i]) {
1100 j = i;
1101 break;
1102 }
1103 }
1104 if (j == -1)
1105 j = narray;
1106 if (j - 1 == lll)
1107 return;
1108
1109 if (j - lll < 0) {
1110 npt = lll - j;
1111 for (i = 0; i < npt; i++) {
1112 array[lll-i] = array[lll-i-1];
1113 iarray[lll-i] = iarray[lll-i-1];
1114 iarray2[lll-i] = iarray2[lll-i-1];
1115 }
1116 array[j] = test;
1117 iarray[j] = it;
1118 iarray2[j] = it2;
1119 } else {
1120 j--;
1121 npt = j - lll;
1122 if (npt != 0) {
1123 for (i = 0; i < npt; i++) {
1124 array[lll+i] = array[lll+i+1];
1125 iarray[lll+i] = iarray[lll+i+1];
1126 iarray2[lll+i] = iarray2[lll+i+1];
1127 }
1128 }
1129 array[j] = test;
1130 iarray[j] = it;
1131 iarray2[j] = it2;
1132 }
1133}
1134
1135/*---------------------------------------------------------------------------*/
1157/*---------------------------------------------------------------------------*/
1158
1159static void sortm(float *a1, unsigned char *a2, int *a3, int n) {
1160 int iii,ii,i,ifin,j,b3;
1161 unsigned char b2;
1162 float b1;
1163
1164 iii = 4;
1165 while (iii < n)
1166 iii *= 2;
1167 iii = min(n,(3*iii)/4 - 1);
1168
1169 while (iii > 1) {
1170 iii /= 2;
1171 ifin = n - iii;
1172 for (ii = 0; ii < ifin; ii++) {
1173 i = ii;
1174 j = i + iii;
1175 if (a1[i] > a1[j]) {
1176 b1 = a1[j];
1177 b2 = a2[j];
1178 b3 = a3[j];
1179 while (1) {
1180 a1[j] = a1[i];
1181 a2[j] = a2[i];
1182 a3[j] = a3[i];
1183 j = i;
1184 i = i - iii;
1185 if (i < 0 || a1[i] <= b1)
1186 break;
1187 }
1188 a1[j] = b1;
1189 a2[j] = b2;
1190 a3[j] = b3;
1191 }
1192 }
1193 }
1194}
1195
1196/*---------------------------------------------------------------------------*/
1216/*---------------------------------------------------------------------------*/
1217
1218static void plugholes(float *data, unsigned char *bpm, int nx) {
1219 int i,ifirst,ilast,i1,i2,j;
1220 float nc,d1,d2,t1,t2,slope;
1221
1222 /* First of all, find the first good value in the array */
1223
1224 i = 0;
1225 while (i < nx && bpm[i] != 0)
1226 i++;
1227 ifirst = i;
1228
1229 /* If all the values in the array are bad, then do nothing */
1230
1231 if (ifirst == nx)
1232 return;
1233
1234 /* Find the last good value in the array */
1235
1236 i = nx - 1;
1237 while (i >= 0 && bpm[i] != 0)
1238 i--;
1239 ilast = i;
1240
1241 /* Right, now start from the first good value and fill in any holes in the
1242 middle part of the array */
1243
1244 i = ifirst;
1245 while (i <= ilast) {
1246 if (bpm[i] == 0) {
1247 i++;
1248 continue;
1249 }
1250 i1 = i - 1;
1251 while (bpm[i] != 0)
1252 i++;
1253 i2 = i;
1254 nc = (float)(i2 - i1 + 1);
1255 d1 = data[i1];
1256 d2 = data[i2];
1257 for (j = i1+1; j <= i2-1; j++) {
1258 t1 = 1.0 - (float)(j - i1)/nc;
1259 t2 = 1.0 - t1;
1260 data[j] = t1*d1 + t2*d2;
1261 }
1262 }
1263
1264 /* Now the left bit... */
1265
1266 if (ifirst > 0) {
1267 slope = data[ifirst+1] - data[ifirst];
1268 for (j = 0; j < ifirst; j++)
1269 data[j] = slope*(float)(j - ifirst) + data[ifirst];
1270 }
1271
1272 /* Now the right bit... */
1273
1274 if (ilast < nx - 1) {
1275 slope = data[ilast] - data[ilast-1];
1276 for (j = ilast; j < nx; j++)
1277 data[j] = slope*(float)(j - ilast) + data[ilast];
1278 }
1279}
1280
1283/*
1284
1285$Log: casu_nebuliser.c,v $
1286Revision 1.3 2015/11/25 10:26:31 jim
1287replaced some hardcoded numbers with defines
1288
1289Revision 1.2 2015/08/07 13:06:54 jim
1290Fixed copyright to ESO
1291
1292Revision 1.1.1.1 2015/06/12 10:44:32 jim
1293Initial import
1294
1295Revision 1.8 2015/03/12 09:16:51 jim
1296Modified to remove some compiler moans
1297
1298Revision 1.7 2015/01/29 11:51:56 jim
1299modified comments
1300
1301Revision 1.6 2014/12/11 12:23:33 jim
1302new version
1303
1304Revision 1.5 2014/05/13 10:26:52 jim
1305Fixed bug in call to medavg
1306
1307Revision 1.4 2014/03/26 15:54:19 jim
1308Removed globals. Confidence is now floating point
1309
1310Revision 1.3 2013/11/21 09:38:14 jim
1311detabbed
1312
1313Revision 1.2 2013-10-24 09:24:33 jim
1314fixed some docs
1315
1316Revision 1.1 2013-09-30 18:09:43 jim
1317Initial entry
1318
1319
1320*/
cpl_image * casu_fits_get_image(casu_fits *p)
Definition: casu_fits.c:436
casu_fits * casu_fits_duplicate(casu_fits *in)
Definition: casu_fits.c:225
cpl_propertylist * casu_fits_get_ehu(casu_fits *p)
Definition: casu_fits.c:576
int casu_nebuliser(casu_fits *infile, casu_fits *inconf, int medfilt, int linfilt, int niter, int axis, int twod, int takeout_sky, int norm, int wantback, float signeg, float sigpos, casu_fits **backmap, int *status)
Remove small scale background variations.
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
float casu_mean(float *data, unsigned char *bpm, long npts)
Definition: casu_stats.c:479