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