/
overscan.py
1339 lines (1180 loc) · 52.7 KB
/
overscan.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
# This file is part of ip_isr.
#
# Developed for the LSST Data Management System.
# This product includes software developed by the LSST Project
# (https://www.lsst.org).
# See the COPYRIGHT file at the top-level directory of this distribution
# for details of code ownership.
#
# This program is free software: you can redistribute it and/or modify
# it under the terms of the GNU General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# This program is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU General Public License for more details.
#
# You should have received a copy of the GNU General Public License
# along with this program. If not, see <https://www.gnu.org/licenses/>.
__all__ = ["OverscanCorrectionTaskConfig", "OverscanCorrectionTask",
"SerialOverscanCorrectionTaskConfig", "SerialOverscanCorrectionTask",
"ParallelOverscanCorrectionTaskConfig", "ParallelOverscanCorrectionTask",
]
import numpy as np
import lsst.afw.math as afwMath
import lsst.afw.image as afwImage
import lsst.geom as geom
import lsst.pipe.base as pipeBase
import lsst.pex.config as pexConfig
from .isr import fitOverscanImage
from .isrFunctions import makeThresholdMask
class OverscanCorrectionTaskConfigBase(pexConfig.Config):
"""Overscan correction options.
"""
fitType = pexConfig.ChoiceField(
dtype=str,
doc="The method for fitting the overscan bias level.",
default='MEDIAN',
allowed={
"POLY": "Fit ordinary polynomial to the longest axis of the overscan region",
"CHEB": "Fit Chebyshev polynomial to the longest axis of the overscan region",
"LEG": "Fit Legendre polynomial to the longest axis of the overscan region",
"NATURAL_SPLINE": "Fit natural spline to the longest axis of the overscan region",
"CUBIC_SPLINE": "Fit cubic spline to the longest axis of the overscan region",
"AKIMA_SPLINE": "Fit Akima spline to the longest axis of the overscan region",
"MEAN": "Correct using the mean of the overscan region",
"MEANCLIP": "Correct using a clipped mean of the overscan region",
"MEDIAN": "Correct using the median of the overscan region",
"MEDIAN_PER_ROW": "Correct using the median per row of the overscan region",
},
)
order = pexConfig.Field(
dtype=int,
doc=("Order of polynomial to fit if overscan fit type is a polynomial, "
"or number of spline knots if overscan fit type is a spline."),
default=1,
)
numSigmaClip = pexConfig.Field(
dtype=float,
doc="Rejection threshold (sigma) for collapsing overscan before fit",
default=3.0,
)
maskPlanes = pexConfig.ListField(
dtype=str,
doc="Mask planes to reject when measuring overscan",
default=['BAD', 'SAT'],
)
overscanIsInt = pexConfig.Field(
dtype=bool,
doc="Treat overscan as an integer image for purposes of fitType=MEDIAN"
" and fitType=MEDIAN_PER_ROW.",
default=True,
)
maxDeviation = pexConfig.Field(
dtype=float,
doc="Maximum deviation from median (in ADU) to mask in overscan correction.",
default=1000.0, check=lambda x: x > 0,
)
class OverscanCorrectionTaskConfig(OverscanCorrectionTaskConfigBase):
doParallelOverscan = pexConfig.Field(
dtype=bool,
doc="Correct using parallel overscan after serial overscan correction?",
default=False,
)
parallelOverscanMaskThreshold = pexConfig.Field(
dtype=int,
doc="Threshold above which pixels in the parallel overscan are masked as bleeds.",
default=100000,
)
parallelOverscanMaskGrowSize = pexConfig.Field(
dtype=int,
doc="Masks created from saturated bleeds should be grown by this many "
"pixels during construction of the parallel overscan mask. "
"This value determined from the ITL chip in the LATISS camera",
default=7,
)
leadingColumnsToSkip = pexConfig.Field(
dtype=int,
doc="Number of leading columns to skip in serial overscan correction.",
default=0,
)
trailingColumnsToSkip = pexConfig.Field(
dtype=int,
doc="Number of trailing columns to skip in serial overscan correction.",
default=0,
)
leadingRowsToSkip = pexConfig.Field(
dtype=int,
doc="Number of leading rows to skip in parallel overscan correction.",
default=0,
)
trailingRowsToSkip = pexConfig.Field(
dtype=int,
doc="Number of trailing rows to skip in parallel overscan correction.",
default=0,
)
class OverscanCorrectionTaskBase(pipeBase.Task):
"""Base Correction task for overscan.
This class contains a number of utilities that are easier to
understand and use when they are not embedded in nested if/else
loops.
Parameters
----------
statControl : `lsst.afw.math.StatisticsControl`, optional
Statistics control object.
"""
ConfigClass = OverscanCorrectionTaskConfigBase
_DefaultName = "overscanBase"
def __init__(self, statControl=None, **kwargs):
super().__init__(**kwargs)
self.allowDebug = True
if statControl:
self.statControl = statControl
else:
self.statControl = afwMath.StatisticsControl()
self.statControl.setNumSigmaClip(self.config.numSigmaClip)
self.statControl.setAndMask(afwImage.Mask.getPlaneBitMask(self.config.maskPlanes))
def run(self, exposure, amp, isTransposed=False):
raise NotImplementedError("run method is not defined for OverscanCorrectionTaskBase")
def correctOverscan(self, exposure, amp, imageBBox, overscanBBox,
isTransposed=True, leadingToSkip=0, trailingToSkip=0):
"""Trim the exposure, fit the overscan, subtract the fit, and
calculate statistics.
Parameters
----------
exposure : `lsst.afw.image.Exposure`
Exposure containing the data.
amp : `lsst.afw.cameraGeom.Amplifier`
The amplifier that is to be corrected.
imageBBox: `lsst.geom.Box2I`
Bounding box of the image data that will have the overscan
subtracted. If parallel overscan will be performed, that
area is added to the image bounding box during serial
overscan correction.
overscanBBox: `lsst.geom.Box2I`
Bounding box for the overscan data.
isTransposed: `bool`
If true, then the data will be transposed before fitting
the overscan.
leadingToSkip : `int`, optional
Leading rows/columns to skip.
trailingToSkip : `int`, optional
Leading rows/columns to skip.
Returns
-------
results : `lsst.pipe.base.Struct`
``ampOverscanModel``
Overscan model broadcast to the full image size.
(`lsst.afw.image.Exposure`)
``overscanOverscanModel``
Overscan model broadcast to the full overscan image
size. (`lsst.afw.image.Exposure`)
``overscanImage``
Overscan image with the overscan fit subtracted.
(`lsst.afw.image.Exposure`)
``overscanValue``
Overscan model. (`float` or `np.array`)
``overscanMean``
Mean value of the overscan fit. (`float`)
``overscanMedian``
Median value of the overscan fit. (`float`)
``overscanSigma``
Standard deviation of the overscan fit. (`float`)
``overscanMeanResidual``
Mean value of the overscan region after overscan
subtraction. (`float`)
``overscanMedianResidual``
Median value of the overscan region after overscan
subtraction. (`float`)
``overscanSigmaResidual``
Standard deviation of the overscan region after
overscan subtraction. (`float`)
"""
overscanBox = self.trimOverscan(exposure, amp, overscanBBox,
leadingToSkip,
trailingToSkip,
transpose=isTransposed)
overscanImage = exposure[overscanBox].getMaskedImage()
overscanArray = overscanImage.image.array
# Mask pixels.
maskVal = overscanImage.mask.getPlaneBitMask(self.config.maskPlanes)
overscanMask = ~((overscanImage.mask.array & maskVal) == 0)
median = np.ma.median(np.ma.masked_where(overscanMask, overscanArray))
bad = np.where(np.abs(overscanArray - median) > self.config.maxDeviation)
overscanImage.mask.array[bad] = overscanImage.mask.getPlaneBitMask("SAT")
# Do overscan fit.
# CZW: Handle transposed correctly.
overscanResults = self.fitOverscan(overscanImage, isTransposed=isTransposed)
# Correct image region (and possibly parallel-overscan region).
ampImage = exposure[imageBBox]
ampOverscanModel = self.broadcastFitToImage(overscanResults.overscanValue,
ampImage.image.array,
transpose=isTransposed)
ampImage.image.array -= ampOverscanModel
# Correct overscan region (and possibly doubly-overscaned
# region).
overscanImage = exposure[overscanBBox]
# CZW: Transposed?
overscanOverscanModel = self.broadcastFitToImage(overscanResults.overscanValue,
overscanImage.image.array)
self.debugView(overscanImage, overscanResults.overscanValue, amp, isTransposed=isTransposed)
overscanImage.image.array -= overscanOverscanModel
# Find residual fit statistics.
stats = afwMath.makeStatistics(overscanImage.getMaskedImage(),
afwMath.MEAN | afwMath.MEDIAN | afwMath.STDEVCLIP, self.statControl)
residualMean = stats.getValue(afwMath.MEAN)
residualMedian = stats.getValue(afwMath.MEDIAN)
residualSigma = stats.getValue(afwMath.STDEVCLIP)
return pipeBase.Struct(ampOverscanModel=ampOverscanModel,
overscanOverscanModel=overscanOverscanModel,
overscanImage=overscanImage,
overscanValue=overscanResults.overscanValue,
overscanMean=overscanResults.overscanMean,
overscanMedian=overscanResults.overscanMedian,
overscanSigma=overscanResults.overscanSigma,
overscanMeanResidual=residualMean,
overscanMedianResidual=residualMedian,
overscanSigmaResidual=residualSigma
)
def broadcastFitToImage(self, overscanValue, imageArray, transpose=False):
"""Broadcast 0 or 1 dimension fit to appropriate shape.
Parameters
----------
overscanValue : `numpy.ndarray`, (Nrows, ) or scalar
Overscan fit to broadcast.
imageArray : `numpy.ndarray`, (Nrows, Ncols)
Image array that we want to match.
transpose : `bool`, optional
Switch order to broadcast along the other axis.
Returns
-------
overscanModel : `numpy.ndarray`, (Nrows, Ncols) or scalar
Expanded overscan fit.
Raises
------
RuntimeError
Raised if no axis has the appropriate dimension.
"""
if isinstance(overscanValue, np.ndarray):
overscanModel = np.zeros_like(imageArray)
if transpose is False:
if imageArray.shape[0] == overscanValue.shape[0]:
overscanModel[:, :] = overscanValue[:, np.newaxis]
elif imageArray.shape[1] == overscanValue.shape[0]:
overscanModel[:, :] = overscanValue[np.newaxis, :]
elif imageArray.shape[0] == overscanValue.shape[1]:
overscanModel[:, :] = overscanValue[np.newaxis, :]
else:
raise RuntimeError(f"Could not broadcast {overscanValue.shape} to "
f"match {imageArray.shape}")
else:
if imageArray.shape[1] == overscanValue.shape[0]:
overscanModel[:, :] = overscanValue[np.newaxis, :]
elif imageArray.shape[0] == overscanValue.shape[0]:
overscanModel[:, :] = overscanValue[:, np.newaxis]
elif imageArray.shape[1] == overscanValue.shape[1]:
overscanModel[:, :] = overscanValue[:, np.newaxis]
else:
raise RuntimeError(f"Could not broadcast {overscanValue.shape} to "
f"match {imageArray.shape}")
else:
overscanModel = overscanValue
return overscanModel
def trimOverscan(self, exposure, amp, bbox, skipLeading, skipTrailing, transpose=False):
"""Trim overscan region to remove edges.
Parameters
----------
exposure : `lsst.afw.image.Exposure`
Exposure containing data.
amp : `lsst.afw.cameraGeom.Amplifier`
Amplifier containing geometry information.
bbox : `lsst.geom.Box2I`
Bounding box of the overscan region.
skipLeading : `int`
Number of leading (towards data region) rows/columns to skip.
skipTrailing : `int`
Number of trailing (away from data region) rows/columns to skip.
transpose : `bool`, optional
Operate on the transposed array.
Returns
-------
overscanArray : `numpy.array`, (N, M)
Data array to fit.
overscanMask : `numpy.array`, (N, M)
Data mask.
"""
dx0, dy0, dx1, dy1 = (0, 0, 0, 0)
dataBBox = amp.getRawDataBBox()
if transpose:
if dataBBox.getBeginY() < bbox.getBeginY():
dy0 += skipLeading
dy1 -= skipTrailing
else:
dy0 += skipTrailing
dy1 -= skipLeading
else:
if dataBBox.getBeginX() < bbox.getBeginX():
dx0 += skipLeading
dx1 -= skipTrailing
else:
dx0 += skipTrailing
dx1 -= skipLeading
overscanBBox = geom.Box2I(bbox.getBegin() + geom.Extent2I(dx0, dy0),
geom.Extent2I(bbox.getWidth() - dx0 + dx1,
bbox.getHeight() - dy0 + dy1))
return overscanBBox
def fitOverscan(self, overscanImage, isTransposed=False):
if self.config.fitType in ('MEAN', 'MEANCLIP', 'MEDIAN'):
# Transposition has no effect here.
overscanResult = self.measureConstantOverscan(overscanImage)
overscanValue = overscanResult.overscanValue
overscanMean = overscanValue
overscanMedian = overscanValue
overscanSigma = 0.0
elif self.config.fitType in ('MEDIAN_PER_ROW', 'POLY', 'CHEB', 'LEG',
'NATURAL_SPLINE', 'CUBIC_SPLINE', 'AKIMA_SPLINE'):
# Force transposes as needed
overscanResult = self.measureVectorOverscan(overscanImage, isTransposed)
overscanValue = overscanResult.overscanValue
stats = afwMath.makeStatistics(overscanResult.overscanValue,
afwMath.MEAN | afwMath.MEDIAN | afwMath.STDEVCLIP,
self.statControl)
overscanMean = stats.getValue(afwMath.MEAN)
overscanMedian = stats.getValue(afwMath.MEDIAN)
overscanSigma = stats.getValue(afwMath.STDEVCLIP)
else:
raise ValueError('%s : %s an invalid overscan type' %
("overscanCorrection", self.config.fitType))
return pipeBase.Struct(overscanValue=overscanValue,
overscanMean=overscanMean,
overscanMedian=overscanMedian,
overscanSigma=overscanSigma,
)
@staticmethod
def integerConvert(image):
"""Return an integer version of the input image.
Parameters
----------
image : `numpy.ndarray`, `lsst.afw.image.Image` or `MaskedImage`
Image to convert to integers.
Returns
-------
outI : `numpy.ndarray`, `lsst.afw.image.Image` or `MaskedImage`
The integer converted image.
Raises
------
RuntimeError
Raised if the input image could not be converted.
"""
if hasattr(image, "image"):
# Is a maskedImage:
imageI = image.image.convertI()
outI = afwImage.MaskedImageI(imageI, image.mask, image.variance)
elif hasattr(image, "convertI"):
# Is an Image:
outI = image.convertI()
elif hasattr(image, "astype"):
# Is a numpy array:
outI = image.astype(int)
else:
raise RuntimeError("Could not convert this to integers: %s %s %s",
image, type(image), dir(image))
return outI
def maskParallelOverscan(self, exposure, detector):
"""Mask the union of high values on all amplifiers in the parallel
overscan.
This operates on the image in-place.
Parameters
----------
exposure : `lsst.afw.image.Exposure`
An untrimmed raw exposure.
detector : `lsst.afw.cameraGeom.Detector`
The detetor to use for amplifier geometry.
"""
parallelMask = None
for amp in detector:
dataView = afwImage.MaskedImageF(exposure.getMaskedImage(),
amp.getRawParallelOverscanBBox(),
afwImage.PARENT)
makeThresholdMask(
maskedImage=dataView,
threshold=self.config.parallelOverscanMaskThreshold,
growFootprints=self.config.parallelOverscanMaskGrowSize,
maskName="BAD"
)
if parallelMask is None:
parallelMask = dataView.mask.array
else:
parallelMask |= dataView.mask.array
for amp in detector:
dataView = afwImage.MaskedImageF(exposure.getMaskedImage(),
amp.getRawParallelOverscanBBox(),
afwImage.PARENT)
dataView.mask.array |= parallelMask
# Constant methods
def measureConstantOverscan(self, image):
"""Measure a constant overscan value.
Parameters
----------
image : `lsst.afw.image.Image` or `lsst.afw.image.MaskedImage`
Image data to measure the overscan from.
Returns
-------
results : `lsst.pipe.base.Struct`
Overscan result with entries:
- ``overscanValue``: Overscan value to subtract (`float`)
- ``isTransposed``: Orientation of the overscan (`bool`)
"""
fitType = afwMath.stringToStatisticsProperty(self.config.fitType)
overscanValue = afwMath.makeStatistics(image, fitType, self.statControl).getValue()
return pipeBase.Struct(overscanValue=overscanValue,
isTransposed=False)
# Vector correction utilities
def getImageArray(self, image):
"""Extract the numpy array from the input image.
Parameters
----------
image : `lsst.afw.image.Image` or `lsst.afw.image.MaskedImage`
Image data to pull array from.
calcImage : `numpy.ndarray`
Image data array for numpy operating.
"""
if hasattr(image, "getImage"):
calcImage = image.getImage().getArray()
calcImage = np.ma.masked_where(image.getMask().getArray() & self.statControl.getAndMask(),
calcImage)
else:
calcImage = image.getArray()
return calcImage
def maskOutliers(self, imageArray):
"""Mask outliers in a row of overscan data from a robust sigma
clipping procedure.
Parameters
----------
imageArray : `numpy.ndarray`
Image to filter along numpy axis=1.
Returns
-------
maskedArray : `numpy.ma.masked_array`
Masked image marking outliers.
"""
lq, median, uq = np.percentile(np.ma.getdata(imageArray),
[25.0, 50.0, 75.0], axis=1)
axisMedians = median
axisStdev = 0.74*(uq - lq) # robust stdev
# Replace pixels that have excessively large stdev values
# with the median of stdev values. A large stdev likely
# indicates a bleed is spilling into the overscan.
axisStdev = np.where(axisStdev > 2.0 * np.median(axisStdev),
np.median(axisStdev), axisStdev)
# Mask pixels that are N-sigma away from their array medians.
diff = np.abs(imageArray - axisMedians[:, np.newaxis])
masked = np.ma.masked_where(diff > self.statControl.getNumSigmaClip()
* axisStdev[:, np.newaxis], imageArray)
return masked
def fillMaskedPixels(self, overscanVector):
"""Fill masked/NaN pixels in the overscan.
Parameters
----------
overscanVector : `np.array` or `np.ma.masked_array`
Overscan vector to fill.
Returns
-------
overscanVector : `np.ma.masked_array`
Filled vector.
Notes
-----
Each maskSlice is a section of overscan with contiguous masks.
Ideally this adds 5 pixels from the left and right of that
mask slice, and takes the median of those values to fill the
slice. If this isn't possible, the median of all non-masked
values is used. The mask is removed for the pixels filled.
"""
workingCopy = overscanVector
if not isinstance(overscanVector, np.ma.MaskedArray):
workingCopy = np.ma.masked_array(overscanVector,
mask=~np.isfinite(overscanVector))
defaultValue = np.median(workingCopy.data[~workingCopy.mask])
for maskSlice in np.ma.clump_masked(workingCopy):
neighborhood = []
if maskSlice.start > 5:
neighborhood.extend(workingCopy[maskSlice.start - 5:maskSlice.start].data)
if maskSlice.stop < workingCopy.size - 5:
neighborhood.extend(workingCopy[maskSlice.stop:maskSlice.stop+5].data)
if len(neighborhood) > 0:
workingCopy.data[maskSlice] = np.nanmedian(neighborhood)
workingCopy.mask[maskSlice] = False
else:
workingCopy.data[maskSlice] = defaultValue
workingCopy.mask[maskSlice] = False
return workingCopy
def collapseArray(self, maskedArray, fillMasked=True):
"""Collapse overscan array (and mask) to a 1-D vector of values.
Parameters
----------
maskedArray : `numpy.ma.masked_array`
Masked array of input overscan data.
fillMasked : `bool`, optional
If true, fill any pixels that are masked with a median of
neighbors.
Returns
-------
collapsed : `numpy.ma.masked_array`
Single dimensional overscan data, combined with the mean.
"""
collapsed = np.mean(maskedArray, axis=1)
if collapsed.mask.sum() > 0 and fillMasked:
collapsed = self.fillMaskedPixels(collapsed)
return collapsed
def collapseArrayMedian(self, maskedArray):
"""Collapse overscan array (and mask) to a 1-D vector of using the
correct integer median of row-values.
Parameters
----------
maskedArray : `numpy.ma.masked_array`
Masked array of input overscan data.
Returns
-------
collapsed : `numpy.ma.masked_array`
Single dimensional overscan data, combined with the afwMath median.
"""
integerMI = self.integerConvert(maskedArray)
collapsed = []
fitType = afwMath.stringToStatisticsProperty('MEDIAN')
for row in integerMI:
newRow = row.compressed()
if len(newRow) > 0:
rowMedian = afwMath.makeStatistics(newRow, fitType, self.statControl).getValue()
else:
rowMedian = np.nan
collapsed.append(rowMedian)
return np.array(collapsed)
def splineFit(self, indices, collapsed, numBins):
"""Wrapper function to match spline fit API to polynomial fit API.
Parameters
----------
indices : `numpy.ndarray`
Locations to evaluate the spline.
collapsed : `numpy.ndarray`
Collapsed overscan values corresponding to the spline
evaluation points.
numBins : `int`
Number of bins to use in constructing the spline.
Returns
-------
interp : `lsst.afw.math.Interpolate`
Interpolation object for later evaluation.
"""
if not np.ma.is_masked(collapsed):
collapsed.mask = np.array(len(collapsed)*[np.ma.nomask])
numPerBin, binEdges = np.histogram(indices, bins=numBins,
weights=1 - collapsed.mask.astype(int))
with np.errstate(invalid="ignore"):
values = np.histogram(indices, bins=numBins,
weights=collapsed.data*~collapsed.mask)[0]/numPerBin
binCenters = np.histogram(indices, bins=numBins,
weights=indices*~collapsed.mask)[0]/numPerBin
if len(binCenters[numPerBin > 0]) < 5:
self.log.warning("Cannot do spline fitting for overscan: %s valid points.",
len(binCenters[numPerBin > 0]))
# Return a scalar value if we have one, otherwise
# return zero. This amplifier is hopefully already
# masked.
if len(values[numPerBin > 0]) != 0:
return float(values[numPerBin > 0][0])
else:
return 0.0
interp = afwMath.makeInterpolate(binCenters.astype(float)[numPerBin > 0],
values.astype(float)[numPerBin > 0],
afwMath.stringToInterpStyle(self.config.fitType))
return interp
@staticmethod
def splineEval(indices, interp):
"""Wrapper function to match spline evaluation API to polynomial fit
API.
Parameters
----------
indices : `numpy.ndarray`
Locations to evaluate the spline.
interp : `lsst.afw.math.interpolate`
Interpolation object to use.
Returns
-------
values : `numpy.ndarray`
Evaluated spline values at each index.
"""
return interp.interpolate(indices.astype(float))
@staticmethod
def maskExtrapolated(collapsed):
"""Create mask if edges are extrapolated.
Parameters
----------
collapsed : `numpy.ma.masked_array`
Masked array to check the edges of.
Returns
-------
maskArray : `numpy.ndarray`
Boolean numpy array of pixels to mask.
"""
maskArray = np.full_like(collapsed, False, dtype=bool)
if np.ma.is_masked(collapsed):
num = len(collapsed)
for low in range(num):
if not collapsed.mask[low]:
break
if low > 0:
maskArray[:low] = True
for high in range(1, num):
if not collapsed.mask[-high]:
break
if high > 1:
maskArray[-high:] = True
return maskArray
def measureVectorOverscan(self, image, isTransposed=False):
"""Calculate the 1-d vector overscan from the input overscan image.
Parameters
----------
image : `lsst.afw.image.MaskedImage`
Image containing the overscan data.
isTransposed : `bool`
If true, the image has been transposed.
Returns
-------
results : `lsst.pipe.base.Struct`
Overscan result with entries:
``overscanValue``
Overscan value to subtract (`float`)
``maskArray``
List of rows that should be masked as ``SUSPECT`` when the
overscan solution is applied. (`list` [ `bool` ])
``isTransposed``
Indicates if the overscan data was transposed during
calcuation, noting along which axis the overscan should be
subtracted. (`bool`)
"""
calcImage = self.getImageArray(image)
# operate on numpy-arrays from here
if isTransposed:
calcImage = np.transpose(calcImage)
masked = self.maskOutliers(calcImage)
if self.config.fitType == 'MEDIAN_PER_ROW':
mi = afwImage.MaskedImageI(image.getBBox())
masked = masked.astype(int)
if isTransposed:
masked = masked.transpose()
mi.image.array[:, :] = masked.data[:, :]
if bool(masked.mask.shape):
mi.mask.array[:, :] = masked.mask[:, :]
overscanVector = fitOverscanImage(mi, self.config.maskPlanes, isTransposed)
overscanVector = self.fillMaskedPixels(overscanVector)
maskArray = self.maskExtrapolated(overscanVector)
else:
collapsed = self.collapseArray(masked)
num = len(collapsed)
indices = 2.0*np.arange(num)/float(num) - 1.0
poly = np.polynomial
fitter, evaler = {
'POLY': (poly.polynomial.polyfit, poly.polynomial.polyval),
'CHEB': (poly.chebyshev.chebfit, poly.chebyshev.chebval),
'LEG': (poly.legendre.legfit, poly.legendre.legval),
'NATURAL_SPLINE': (self.splineFit, self.splineEval),
'CUBIC_SPLINE': (self.splineFit, self.splineEval),
'AKIMA_SPLINE': (self.splineFit, self.splineEval)
}[self.config.fitType]
# These are the polynomial coefficients, or an
# interpolation object.
coeffs = fitter(indices, collapsed, self.config.order)
if isinstance(coeffs, float):
self.log.warning("Using fallback value %f due to fitter failure. Amplifier will be masked.",
coeffs)
overscanVector = np.full_like(indices, coeffs)
maskArray = np.full_like(collapsed, True, dtype=bool)
else:
# Otherwise we can just use things as normal.
overscanVector = evaler(indices, coeffs)
maskArray = self.maskExtrapolated(collapsed)
return pipeBase.Struct(overscanValue=np.array(overscanVector),
maskArray=maskArray,
isTransposed=isTransposed)
def debugView(self, image, model, amp=None, isTransposed=True):
"""Debug display for the final overscan solution.
Parameters
----------
image : `lsst.afw.image.Image`
Input image the overscan solution was determined from.
model : `numpy.ndarray` or `float`
Overscan model determined for the image.
amp : `lsst.afw.cameraGeom.Amplifier`, optional
Amplifier to extract diagnostic information.
isTransposed : `bool`, optional
Does the data need to be transposed before display?
"""
import lsstDebug
if not lsstDebug.Info(__name__).display:
return
if not self.allowDebug:
return
calcImage = self.getImageArray(image)
# CZW: Check that this is ok
if isTransposed:
calcImage = np.transpose(calcImage)
masked = self.maskOutliers(calcImage)
collapsed = self.collapseArray(masked, fillMasked=False)
num = len(collapsed)
indices = 2.0 * np.arange(num)/float(num) - 1.0
indices = np.arange(num)
if np.ma.is_masked(collapsed):
collapsedMask = collapsed.mask
else:
collapsedMask = np.array(num*[np.ma.nomask])
import matplotlib.pyplot as plot
figure = plot.figure(1)
figure.clear()
axes = figure.add_axes((0.1, 0.1, 0.8, 0.8))
axes.plot(indices[~collapsedMask], collapsed[~collapsedMask], 'k+')
if collapsedMask.sum() > 0:
axes.plot(indices[collapsedMask], collapsed.data[collapsedMask], 'b+')
if isinstance(model, np.ndarray):
plotModel = model
else:
plotModel = np.zeros_like(indices)
plotModel += model
axes.plot(indices, plotModel, 'r-')
plot.xlabel("position along overscan region")
plot.ylabel("pixel value/fit value")
if amp:
plot.title(f"{amp.getName()} DataX: "
f"[{amp.getRawDataBBox().getBeginX()}:{amp.getRawBBox().getEndX()}]"
f"OscanX: [{amp.getRawHorizontalOverscanBBox().getBeginX()}:"
f"{amp.getRawHorizontalOverscanBBox().getEndX()}] {self.config.fitType}")
else:
plot.title("No amp supplied.")
figure.show()
prompt = "Press Enter or c to continue [chp]..."
while True:
ans = input(prompt).lower()
if ans in ("", " ", "c",):
break
elif ans in ("p", ):
import pdb
pdb.set_trace()
elif ans in ('x', ):
self.allowDebug = False
break
elif ans in ("h", ):
print("[h]elp [c]ontinue [p]db e[x]itDebug")
plot.close()
class OverscanCorrectionTask(OverscanCorrectionTaskBase):
"""Correction task for serial/parallel overscan.
(Will be deprecated)
This class contains a number of utilities that are easier to
understand and use when they are not embedded in nested if/else
loops.
Parameters
----------
statControl : `lsst.afw.math.StatisticsControl`, optional
Statistics control object.
"""
ConfigClass = OverscanCorrectionTaskConfig
_DefaultName = "overscan"
def run(self, exposure, amp, isTransposed=False):
"""Measure and remove serial/parallel overscan from an amplifier image.
This will be deprecated.
Parameters
----------
exposure : `lsst.afw.image.Exposure`
Image data that will have the overscan corrections applied.
amp : `lsst.afw.cameraGeom.Amplifier`
Amplifier to use for debugging purposes.
isTransposed : `bool`, optional
Is the image transposed, such that serial and parallel
overscan regions are reversed? Default is False.
Returns
-------
overscanResults : `lsst.pipe.base.Struct`
Result struct with components:
``imageFit``
Value or fit subtracted from the amplifier image data
(scalar or `lsst.afw.image.Image`).
``overscanFit``
Value or fit subtracted from the serial overscan image
data (scalar or `lsst.afw.image.Image`).
``overscanImage``
Image of the serial overscan region with the serial
overscan correction applied
(`lsst.afw.image.Image`). This quantity is used to
estimate the amplifier read noise empirically.
``parallelOverscanFit``
Value or fit subtracted from the parallel overscan
image data (scalar, `lsst.afw.image.Image`, or None).
``parallelOverscanImage``
Image of the parallel overscan region with the
parallel overscan correction applied
(`lsst.afw.image.Image` or None).
``overscanMean``
Mean of the fit serial overscan region.
This and the following values will be tuples of
(serial, parallel) if doParallelOverscan=True.
``overscanMedian``
Median of the fit serial overscan region.
``overscanSigma``
Sigma of the fit serial overscan region.
``residualMean``
Mean of the residual of the serial overscan region after
correction.
``residualMedian``
Median of the residual of the serial overscan region after
correction.
``residualSigma``
Mean of the residual of the serial overscan region after
correction.
Raises
------
RuntimeError
Raised if an invalid overscan type is set.
"""
# Do Serial overscan first.
serialOverscanBBox = amp.getRawSerialOverscanBBox()
imageBBox = amp.getRawDataBBox()
if self.config.doParallelOverscan:
# We need to extend the serial overscan BBox to the full
# size of the detector.
parallelOverscanBBox = amp.getRawParallelOverscanBBox()
imageBBox = imageBBox.expandedTo(parallelOverscanBBox)
if isTransposed:
serialOverscanBBox = geom.Box2I(
geom.Point2I(serialOverscanBBox.getMinX(), imageBBox.getEndY()),
geom.Extent2I(imageBBox.getWidth(), serialOverscanBBox.getHeight()),
)
else:
serialOverscanBBox = geom.Box2I(
geom.Point2I(serialOverscanBBox.getMinX(),
imageBBox.getMinY()),
geom.Extent2I(serialOverscanBBox.getWidth(),
imageBBox.getHeight()),
)
serialResults = self.correctOverscan(
exposure,
amp,
imageBBox,
serialOverscanBBox,
isTransposed=isTransposed,
leadingToSkip=self.config.leadingColumnsToSkip,
trailingToSkip=self.config.trailingColumnsToSkip,
)
overscanMean = serialResults.overscanMean
overscanMedian = serialResults.overscanMedian
overscanSigma = serialResults.overscanSigma
residualMean = serialResults.overscanMeanResidual
residualMedian = serialResults.overscanMedianResidual
residualSigma = serialResults.overscanSigmaResidual
# Do Parallel Overscan
parallelResults = None
if self.config.doParallelOverscan:
# This does not need any extensions, as we'll only
# subtract it from the data region.