File: kernels.py

package info (click to toggle)
python-astropy 1.3-8~bpo8%2B2
  • links: PTS, VCS
  • area: main
  • in suites: jessie-backports
  • size: 44,292 kB
  • sloc: ansic: 160,360; python: 137,322; sh: 11,493; lex: 7,638; yacc: 4,956; xml: 1,796; makefile: 474; cpp: 364
file content (1018 lines) | stat: -rw-r--r-- 32,465 bytes parent folder | download | duplicates (2)
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
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
# Licensed under a 3-clause BSD style license - see LICENSE.rst
from __future__ import (absolute_import, division, print_function,
                        unicode_literals)

import math

import numpy as np

from .core import Kernel1D, Kernel2D, Kernel
from .utils import KernelSizeError
from ..modeling import models
from ..modeling.core import Fittable1DModel, Fittable2DModel


__all__ = ['Gaussian1DKernel', 'Gaussian2DKernel', 'CustomKernel',
           'Box1DKernel', 'Box2DKernel', 'Tophat2DKernel',
           'Trapezoid1DKernel', 'MexicanHat1DKernel', 'MexicanHat2DKernel',
           'AiryDisk2DKernel', 'Moffat2DKernel', 'Model1DKernel',
           'Model2DKernel', 'TrapezoidDisk2DKernel', 'Ring2DKernel']


def _round_up_to_odd_integer(value):
    i = int(math.ceil(value))  # TODO: int() call is only needed for six.PY2
    if i % 2 == 0:
        return i + 1
    else:
        return i


class Gaussian1DKernel(Kernel1D):
    """
    1D Gaussian filter kernel.

    The Gaussian filter is a filter with great smoothing properties. It is
    isotropic and does not produce artifacts.

    Parameters
    ----------
    stddev : number
        Standard deviation of the Gaussian kernel.
    x_size : odd int, optional
        Size of the kernel array. Default = 8 * stddev
    mode : str, optional
        One of the following discretization modes:
            * 'center' (default)
                Discretize model by taking the value
                at the center of the bin.
            * 'linear_interp'
                Discretize model by linearly interpolating
                between the values at the corners of the bin.
            * 'oversample'
                Discretize model by taking the average
                on an oversampled grid.
            * 'integrate'
                Discretize model by integrating the
                model over the bin. Very slow.
    factor : number, optional
        Factor of oversampling. Default factor = 10. If the factor
        is too large, evaluation can be very slow.


    See Also
    --------
    Box1DKernel, Trapezoid1DKernel, MexicanHat1DKernel


    Examples
    --------
    Kernel response:

     .. plot::
        :include-source:

        import matplotlib.pyplot as plt
        from astropy.convolution import Gaussian1DKernel
        gauss_1D_kernel = Gaussian1DKernel(10)
        plt.plot(gauss_1D_kernel, drawstyle='steps')
        plt.xlabel('x [pixels]')
        plt.ylabel('value')
        plt.show()
    """
    _separable = True
    _is_bool = False

    def __init__(self, stddev, **kwargs):
        self._model = models.Gaussian1D(1. / (np.sqrt(2 * np.pi) * stddev),
                                        0, stddev)
        self._default_size = _round_up_to_odd_integer(8 * stddev)
        super(Gaussian1DKernel, self).__init__(**kwargs)
        self._truncation = np.abs(1. - self._array.sum())


class Gaussian2DKernel(Kernel2D):
    """
    2D Gaussian filter kernel.

    The Gaussian filter is a filter with great smoothing properties. It is
    isotropic and does not produce artifacts.

    Parameters
    ----------
    stddev : number
        Standard deviation of the Gaussian kernel.
    x_size : odd int, optional
        Size in x direction of the kernel array. Default = 8 * stddev.
    y_size : odd int, optional
        Size in y direction of the kernel array. Default = 8 * stddev.
    mode : str, optional
        One of the following discretization modes:
            * 'center' (default)
                Discretize model by taking the value
                at the center of the bin.
            * 'linear_interp'
                Discretize model by performing a bilinear interpolation
                between the values at the corners of the bin.
            * 'oversample'
                Discretize model by taking the average
                on an oversampled grid.
            * 'integrate'
                Discretize model by integrating the
                model over the bin.
    factor : number, optional
        Factor of oversampling. Default factor = 10.


    See Also
    --------
    Box2DKernel, Tophat2DKernel, MexicanHat2DKernel, Ring2DKernel,
    TrapezoidDisk2DKernel, AiryDisk2DKernel, Moffat2DKernel

    Examples
    --------
    Kernel response:

     .. plot::
        :include-source:

        import matplotlib.pyplot as plt
        from astropy.convolution import Gaussian2DKernel
        gaussian_2D_kernel = Gaussian2DKernel(10)
        plt.imshow(gaussian_2D_kernel, interpolation='none', origin='lower')
        plt.xlabel('x [pixels]')
        plt.ylabel('y [pixels]')
        plt.colorbar()
        plt.show()

    """
    _separable = True
    _is_bool = False

    def __init__(self, stddev, **kwargs):
        self._model = models.Gaussian2D(1. / (2 * np.pi * stddev ** 2), 0,
                                        0, stddev, stddev)
        self._default_size = _round_up_to_odd_integer(8 * stddev)
        super(Gaussian2DKernel, self).__init__(**kwargs)
        self._truncation = np.abs(1. - self._array.sum())


class Box1DKernel(Kernel1D):
    """
    1D Box filter kernel.

    The Box filter or running mean is a smoothing filter. It is not isotropic
    and can produce artifacts, when applied repeatedly to the same data.

    By default the Box kernel uses the ``linear_interp`` discretization mode,
    which allows non-shifting, even-sized kernels.  This is achieved by
    weighting the edge pixels with 1/2. E.g a Box kernel with an effective
    smoothing of 4 pixel would have the following array: [0.5, 1, 1, 1, 0.5].


    Parameters
    ----------
    width : number
        Width of the filter kernel.
    mode : str, optional
        One of the following discretization modes:
            * 'center'
                Discretize model by taking the value
                at the center of the bin.
            * 'linear_interp' (default)
                Discretize model by linearly interpolating
                between the values at the corners of the bin.
            * 'oversample'
                Discretize model by taking the average
                on an oversampled grid.
            * 'integrate'
                Discretize model by integrating the
                model over the bin.
    factor : number, optional
        Factor of oversampling. Default factor = 10.

    See Also
    --------
    Gaussian1DKernel, Trapezoid1DKernel, MexicanHat1DKernel


    Examples
    --------
    Kernel response function:

     .. plot::
        :include-source:

        import matplotlib.pyplot as plt
        from astropy.convolution import Box1DKernel
        box_1D_kernel = Box1DKernel(9)
        plt.plot(box_1D_kernel, drawstyle='steps')
        plt.xlim(-1, 9)
        plt.xlabel('x [pixels]')
        plt.ylabel('value')
        plt.show()

    """
    _separable = True
    _is_bool = True

    def __init__(self, width, **kwargs):
        self._model = models.Box1D(1. / width, 0, width)
        self._default_size = _round_up_to_odd_integer(width)
        kwargs['mode'] = 'linear_interp'
        super(Box1DKernel, self).__init__(**kwargs)
        self._truncation = 0
        self.normalize()


class Box2DKernel(Kernel2D):
    """
    2D Box filter kernel.

    The Box filter or running mean is a smoothing filter. It is not isotropic
    and can produce artifact, when applied repeatedly to the same data.

    By default the Box kernel uses the ``linear_interp`` discretization mode,
    which allows non-shifting, even-sized kernels.  This is achieved by
    weighting the edge pixels with 1/2.


    Parameters
    ----------
    width : number
        Width of the filter kernel.
    mode : str, optional
        One of the following discretization modes:
            * 'center'
                Discretize model by taking the value
                at the center of the bin.
            * 'linear_interp' (default)
                Discretize model by performing a bilinear interpolation
                between the values at the corners of the bin.
            * 'oversample'
                Discretize model by taking the average
                on an oversampled grid.
            * 'integrate'
                Discretize model by integrating the
                model over the bin.
    factor : number, optional
        Factor of oversampling. Default factor = 10.


    See Also
    --------
    Gaussian2DKernel, Tophat2DKernel, MexicanHat2DKernel, Ring2DKernel,
    TrapezoidDisk2DKernel, AiryDisk2DKernel, Moffat2DKernel

    Examples
    --------
    Kernel response:

     .. plot::
        :include-source:

        import matplotlib.pyplot as plt
        from astropy.convolution import Box2DKernel
        box_2D_kernel = Box2DKernel(9)
        plt.imshow(box_2D_kernel, interpolation='none', origin='lower',
                   vmin=0.0, vmax=0.015)
        plt.xlim(-1, 9)
        plt.ylim(-1, 9)
        plt.xlabel('x [pixels]')
        plt.ylabel('y [pixels]')
        plt.colorbar()
        plt.show()
    """
    _separable = True
    _is_bool = True

    def __init__(self, width, **kwargs):
        self._model = models.Box2D(1. / width ** 2, 0, 0, width, width)
        self._default_size = _round_up_to_odd_integer(width)
        kwargs['mode'] = 'linear_interp'
        super(Box2DKernel, self).__init__(**kwargs)
        self._truncation = 0
        self.normalize()


class Tophat2DKernel(Kernel2D):
    """
    2D Tophat filter kernel.

    The Tophat filter is an isotropic smoothing filter. It can produce
    artifacts when applied repeatedly on the same data.

    Parameters
    ----------
    radius : int
        Radius of the filter kernel.
    mode : str, optional
        One of the following discretization modes:
            * 'center' (default)
                Discretize model by taking the value
                at the center of the bin.
            * 'linear_interp'
                Discretize model by performing a bilinear interpolation
                between the values at the corners of the bin.
            * 'oversample'
                Discretize model by taking the average
                on an oversampled grid.
            * 'integrate'
                Discretize model by integrating the
                model over the bin.
    factor : number, optional
        Factor of oversampling. Default factor = 10.


    See Also
    --------
    Gaussian2DKernel, Box2DKernel, MexicanHat2DKernel, Ring2DKernel,
    TrapezoidDisk2DKernel, AiryDisk2DKernel, Moffat2DKernel

    Examples
    --------
    Kernel response:

     .. plot::
        :include-source:

        import matplotlib.pyplot as plt
        from astropy.convolution import Tophat2DKernel
        tophat_2D_kernel = Tophat2DKernel(40)
        plt.imshow(tophat_2D_kernel, interpolation='none', origin='lower')
        plt.xlabel('x [pixels]')
        plt.ylabel('y [pixels]')
        plt.colorbar()
        plt.show()

    """
    def __init__(self, radius, **kwargs):
        self._model = models.Disk2D(1. / (np.pi * radius ** 2), 0, 0, radius)
        self._default_size = _round_up_to_odd_integer(2 * radius)
        super(Tophat2DKernel, self).__init__(**kwargs)
        self._truncation = 0


class Ring2DKernel(Kernel2D):
    """
    2D Ring filter kernel.

    The Ring filter kernel is the difference between two Tophat kernels of
    different width. This kernel is useful for, e.g., background estimation.

    Parameters
    ----------
    radius_in : number
        Inner radius of the ring kernel.
    width : number
        Width of the ring kernel.
    mode : str, optional
        One of the following discretization modes:
            * 'center' (default)
                Discretize model by taking the value
                at the center of the bin.
            * 'linear_interp'
                Discretize model by performing a bilinear interpolation
                between the values at the corners of the bin.
            * 'oversample'
                Discretize model by taking the average
                on an oversampled grid.
            * 'integrate'
                Discretize model by integrating the
                model over the bin.
    factor : number, optional
        Factor of oversampling. Default factor = 10.

    See Also
    --------
    Gaussian2DKernel, Box2DKernel, Tophat2DKernel, MexicanHat2DKernel,
    Ring2DKernel, AiryDisk2DKernel, Moffat2DKernel

    Examples
    --------
    Kernel response:

     .. plot::
        :include-source:

        import matplotlib.pyplot as plt
        from astropy.convolution import Ring2DKernel
        ring_2D_kernel = Ring2DKernel(9, 8)
        plt.imshow(ring_2D_kernel, interpolation='none', origin='lower')
        plt.xlabel('x [pixels]')
        plt.ylabel('y [pixels]')
        plt.colorbar()
        plt.show()
    """
    def __init__(self, radius_in, width, **kwargs):
        radius_out = radius_in + width
        self._model = models.Ring2D(1. / (np.pi * (radius_out ** 2 - radius_in ** 2)),
                                    0, 0, radius_in, width)
        self._default_size = _round_up_to_odd_integer(2 * radius_out)
        super(Ring2DKernel, self).__init__(**kwargs)
        self._truncation = 0


class Trapezoid1DKernel(Kernel1D):
    """
    1D trapezoid kernel.

    Parameters
    ----------
    width : number
        Width of the filter kernel, defined as the width of the constant part,
        before it begins to slope down.
    slope : number
        Slope of the filter kernel's tails
    mode : str, optional
        One of the following discretization modes:
            * 'center' (default)
                Discretize model by taking the value
                at the center of the bin.
            * 'linear_interp'
                Discretize model by linearly interpolating
                between the values at the corners of the bin.
            * 'oversample'
                Discretize model by taking the average
                on an oversampled grid.
            * 'integrate'
                Discretize model by integrating the
                model over the bin.
    factor : number, optional
        Factor of oversampling. Default factor = 10.

    See Also
    --------
    Box1DKernel, Gaussian1DKernel, MexicanHat1DKernel

    Examples
    --------
    Kernel response:

     .. plot::
        :include-source:

        import matplotlib.pyplot as plt
        from astropy.convolution import Trapezoid1DKernel
        trapezoid_1D_kernel = Trapezoid1DKernel(17, slope=0.2)
        plt.plot(trapezoid_1D_kernel, drawstyle='steps')
        plt.xlabel('x [pixels]')
        plt.ylabel('amplitude')
        plt.xlim(-1, 28)
        plt.show()
    """
    _is_bool = False

    def __init__(self, width, slope=1., **kwargs):
        self._model = models.Trapezoid1D(1, 0, width, slope)
        self._default_size = _round_up_to_odd_integer(width + 2. / slope)
        super(Trapezoid1DKernel, self).__init__(**kwargs)
        self._truncation = 0
        self.normalize()


class TrapezoidDisk2DKernel(Kernel2D):
    """
    2D trapezoid kernel.

    Parameters
    ----------
    radius : number
        Width of the filter kernel, defined as the width of the constant part,
        before it begins to slope down.
    slope : number
        Slope of the filter kernel's tails
    mode : str, optional
        One of the following discretization modes:
            * 'center' (default)
                Discretize model by taking the value
                at the center of the bin.
            * 'linear_interp'
                Discretize model by performing a bilinear interpolation
                between the values at the corners of the bin.
            * 'oversample'
                Discretize model by taking the average
                on an oversampled grid.
            * 'integrate'
                Discretize model by integrating the
                model over the bin.
    factor : number, optional
        Factor of oversampling. Default factor = 10.

    See Also
    --------
    Gaussian2DKernel, Box2DKernel, Tophat2DKernel, MexicanHat2DKernel,
    Ring2DKernel, AiryDisk2DKernel, Moffat2DKernel

    Examples
    --------
    Kernel response:

     .. plot::
        :include-source:

        import matplotlib.pyplot as plt
        from astropy.convolution import TrapezoidDisk2DKernel
        trapezoid_2D_kernel = TrapezoidDisk2DKernel(20, slope=0.2)
        plt.imshow(trapezoid_2D_kernel, interpolation='none', origin='lower')
        plt.xlabel('x [pixels]')
        plt.ylabel('y [pixels]')
        plt.colorbar()
        plt.show()

    """
    _is_bool = False

    def __init__(self, radius, slope=1., **kwargs):
        self._model = models.TrapezoidDisk2D(1, 0, 0, radius, slope)
        self._default_size = _round_up_to_odd_integer(2 * radius + 2. / slope)
        super(TrapezoidDisk2DKernel, self).__init__(**kwargs)
        self._truncation = 0
        self.normalize()


class MexicanHat1DKernel(Kernel1D):
    """
    1D Mexican hat filter kernel.

    The Mexican Hat, or inverted Gaussian-Laplace filter, is a
    bandpass filter. It smoothes the data and removes slowly varying
    or constant structures (e.g. Background). It is useful for peak or
    multi-scale detection.

    This kernel is derived from a normalized Gaussian function, by
    computing the second derivative. This results in an amplitude
    at the kernels center of 1. / (sqrt(2 * pi) * width ** 3). The
    normalization is the same as for `scipy.ndimage.gaussian_laplace`,
    except for a minus sign.

    Parameters
    ----------
    width : number
        Width of the filter kernel, defined as the standard deviation
        of the Gaussian function from which it is derived.
    x_size : odd int, optional
        Size in x direction of the kernel array. Default = 8 * width.
    mode : str, optional
        One of the following discretization modes:
            * 'center' (default)
                Discretize model by taking the value
                at the center of the bin.
            * 'linear_interp'
                Discretize model by linearly interpolating
                between the values at the corners of the bin.
            * 'oversample'
                Discretize model by taking the average
                on an oversampled grid.
            * 'integrate'
                Discretize model by integrating the
                model over the bin.
    factor : number, optional
        Factor of oversampling. Default factor = 10.


    See Also
    --------
    Box1DKernel, Gaussian1DKernel, Trapezoid1DKernel

    Examples
    --------
    Kernel response:

     .. plot::
        :include-source:

        import matplotlib.pyplot as plt
        from astropy.convolution import MexicanHat1DKernel
        mexicanhat_1D_kernel = MexicanHat1DKernel(10)
        plt.plot(mexicanhat_1D_kernel, drawstyle='steps')
        plt.xlabel('x [pixels]')
        plt.ylabel('value')
        plt.show()

    """
    _is_bool = True

    def __init__(self, width, **kwargs):
        amplitude = 1.0 / (np.sqrt(2 * np.pi) * width ** 3)
        self._model = models.MexicanHat1D(amplitude, 0, width)
        self._default_size = _round_up_to_odd_integer(8 * width)
        super(MexicanHat1DKernel, self).__init__(**kwargs)
        self._truncation = np.abs(self._array.sum() / self._array.size)


class MexicanHat2DKernel(Kernel2D):
    """
    2D Mexican hat filter kernel.

    The Mexican Hat, or inverted Gaussian-Laplace filter, is a
    bandpass filter. It smoothes the data and removes slowly varying
    or constant structures (e.g. Background). It is useful for peak or
    multi-scale detection.

    This kernel is derived from a normalized Gaussian function, by
    computing the second derivative. This results in an amplitude
    at the kernels center of 1. / (pi * width ** 4). The normalization
    is the same as for `scipy.ndimage.gaussian_laplace`, except
    for a minus sign.

    Parameters
    ----------
    width : number
        Width of the filter kernel, defined as the standard deviation
        of the Gaussian function from which it is derived.
    x_size : odd int, optional
        Size in x direction of the kernel array. Default = 8 * width.
    y_size : odd int, optional
        Size in y direction of the kernel array. Default = 8 * width.
    mode : str, optional
        One of the following discretization modes:
            * 'center' (default)
                Discretize model by taking the value
                at the center of the bin.
            * 'linear_interp'
                Discretize model by performing a bilinear interpolation
                between the values at the corners of the bin.
            * 'oversample'
                Discretize model by taking the average
                on an oversampled grid.
            * 'integrate'
                Discretize model by integrating the
                model over the bin.
    factor : number, optional
        Factor of oversampling. Default factor = 10.


    See Also
    --------
    Gaussian2DKernel, Box2DKernel, Tophat2DKernel, Ring2DKernel,
    TrapezoidDisk2DKernel, AiryDisk2DKernel, Moffat2DKernel

    Examples
    --------
    Kernel response:

     .. plot::
        :include-source:

        import matplotlib.pyplot as plt
        from astropy.convolution import MexicanHat2DKernel
        mexicanhat_2D_kernel = MexicanHat2DKernel(10)
        plt.imshow(mexicanhat_2D_kernel, interpolation='none', origin='lower')
        plt.xlabel('x [pixels]')
        plt.ylabel('y [pixels]')
        plt.colorbar()
        plt.show()
    """
    _is_bool = False

    def __init__(self, width, **kwargs):
        amplitude = 1.0 / (np.pi * width ** 4)
        self._model = models.MexicanHat2D(amplitude, 0, 0, width)
        self._default_size = _round_up_to_odd_integer(8 * width)
        super(MexicanHat2DKernel, self).__init__(**kwargs)
        self._truncation = np.abs(self._array.sum() / self._array.size)


class AiryDisk2DKernel(Kernel2D):
    """
    2D Airy disk kernel.

    This kernel models the diffraction pattern of a circular aperture. This
    kernel is normalized to a peak value of 1.

    Parameters
    ----------
    radius : float
        The radius of the Airy disk kernel (radius of the first zero).
    x_size : odd int, optional
        Size in x direction of the kernel array. Default = 8 * radius.
    y_size : odd int, optional
        Size in y direction of the kernel array. Default = 8 * radius.
    mode : str, optional
        One of the following discretization modes:
            * 'center' (default)
                Discretize model by taking the value
                at the center of the bin.
            * 'linear_interp'
                Discretize model by performing a bilinear interpolation
                between the values at the corners of the bin.
            * 'oversample'
                Discretize model by taking the average
                on an oversampled grid.
            * 'integrate'
                Discretize model by integrating the
                model over the bin.
    factor : number, optional
        Factor of oversampling. Default factor = 10.

    See Also
    --------
    Gaussian2DKernel, Box2DKernel, Tophat2DKernel, MexicanHat2DKernel,
    Ring2DKernel, TrapezoidDisk2DKernel, AiryDisk2DKernel, Moffat2DKernel

    Examples
    --------
    Kernel response:

     .. plot::
        :include-source:

        import matplotlib.pyplot as plt
        from astropy.convolution import AiryDisk2DKernel
        airydisk_2D_kernel = AiryDisk2DKernel(10)
        plt.imshow(airydisk_2D_kernel, interpolation='none', origin='lower')
        plt.xlabel('x [pixels]')
        plt.ylabel('y [pixels]')
        plt.colorbar()
        plt.show()
    """
    _is_bool = False

    def __init__(self, radius, **kwargs):
        self._model = models.AiryDisk2D(1, 0, 0, radius)
        self._default_size = _round_up_to_odd_integer(8 * radius)
        super(AiryDisk2DKernel, self).__init__(**kwargs)
        self.normalize()
        self._truncation = None


class Moffat2DKernel(Kernel2D):
    """
    2D Moffat kernel.

    This kernel is a typical model for a seeing limited PSF.

    Parameters
    ----------
    gamma : float
        Core width of the Moffat model.
    alpha : float
        Power index of the Moffat model.
    x_size : odd int, optional
        Size in x direction of the kernel array. Default = 8 * radius.
    y_size : odd int, optional
        Size in y direction of the kernel array. Default = 8 * radius.
    mode : str, optional
        One of the following discretization modes:
            * 'center' (default)
                Discretize model by taking the value
                at the center of the bin.
            * 'linear_interp'
                Discretize model by performing a bilinear interpolation
                between the values at the corners of the bin.
            * 'oversample'
                Discretize model by taking the average
                on an oversampled grid.
            * 'integrate'
                Discretize model by integrating the
                model over the bin.
    factor : number, optional
        Factor of oversampling. Default factor = 10.

    See Also
    --------
    Gaussian2DKernel, Box2DKernel, Tophat2DKernel, MexicanHat2DKernel,
    Ring2DKernel, TrapezoidDisk2DKernel, AiryDisk2DKernel

    Examples
    --------
    Kernel response:

     .. plot::
        :include-source:

        import matplotlib.pyplot as plt
        from astropy.convolution import Moffat2DKernel
        moffat_2D_kernel = Moffat2DKernel(3, 2)
        plt.imshow(moffat_2D_kernel, interpolation='none', origin='lower')
        plt.xlabel('x [pixels]')
        plt.ylabel('y [pixels]')
        plt.colorbar()
        plt.show()
    """
    _is_bool = False

    def __init__(self, gamma, alpha, **kwargs):
        self._model = models.Moffat2D((gamma - 1.0) / (np.pi * alpha * alpha),
                                      0, 0, gamma, alpha)
        fwhm = 2.0 * alpha * (2.0 ** (1.0 / gamma) - 1.0) ** 0.5
        self._default_size = _round_up_to_odd_integer(4.0 * fwhm)
        super(Moffat2DKernel, self).__init__(**kwargs)
        self.normalize()
        self._truncation = None


class Model1DKernel(Kernel1D):
    """
    Create kernel from 1D model.

    The model has to be centered on x = 0.

    Parameters
    ----------
    model : `~astropy.modeling.Fittable1DModel`
        Kernel response function model
    x_size : odd int, optional
        Size in x direction of the kernel array. Default = 8 * width.
    mode : str, optional
        One of the following discretization modes:
            * 'center' (default)
                Discretize model by taking the value
                at the center of the bin.
            * 'linear_interp'
                Discretize model by linearly interpolating
                between the values at the corners of the bin.
            * 'oversample'
                Discretize model by taking the average
                on an oversampled grid.
            * 'integrate'
                Discretize model by integrating the
                model over the bin.
    factor : number, optional
        Factor of oversampling. Default factor = 10.

    Raises
    ------
    TypeError
        If model is not an instance of `~astropy.modeling.Fittable1DModel`

    See also
    --------
    Model2DKernel : Create kernel from `~astropy.modeling.Fittable2DModel`
    CustomKernel : Create kernel from list or array

    Examples
    --------
    Define a Gaussian1D model:

        >>> from astropy.modeling.models import Gaussian1D
        >>> from astropy.convolution.kernels import Model1DKernel
        >>> gauss = Gaussian1D(1, 0, 2)

    And create a custom one dimensional kernel from it:

        >>> gauss_kernel = Model1DKernel(gauss, x_size=9)

    This kernel can now be used like a usual Astropy kernel.
    """
    _separable = False
    _is_bool = False

    def __init__(self, model, **kwargs):
        if isinstance(model, Fittable1DModel):
            self._model = model
        else:
            raise TypeError("Must be Fittable1DModel")
        super(Model1DKernel, self).__init__(**kwargs)


class Model2DKernel(Kernel2D):
    """
    Create kernel from 2D model.

    The model has to be centered on x = 0 and y = 0.

    Parameters
    ----------
    model : `~astropy.modeling.Fittable2DModel`
        Kernel response function model
    x_size : odd int, optional
        Size in x direction of the kernel array. Default = 8 * width.
    y_size : odd int, optional
        Size in y direction of the kernel array. Default = 8 * width.
    mode : str, optional
        One of the following discretization modes:
            * 'center' (default)
                Discretize model by taking the value
                at the center of the bin.
            * 'linear_interp'
                Discretize model by performing a bilinear interpolation
                between the values at the corners of the bin.
            * 'oversample'
                Discretize model by taking the average
                on an oversampled grid.
            * 'integrate'
                Discretize model by integrating the
                model over the bin.
    factor : number, optional
        Factor of oversampling. Default factor = 10.

    Raises
    ------
    TypeError
        If model is not an instance of `~astropy.modeling.Fittable2DModel`

    See also
    --------
    Model1DKernel : Create kernel from `~astropy.modeling.Fittable1DModel`
    CustomKernel : Create kernel from list or array

    Examples
    --------
    Define a Gaussian2D model:

        >>> from astropy.modeling.models import Gaussian2D
        >>> from astropy.convolution.kernels import Model2DKernel
        >>> gauss = Gaussian2D(1, 0, 0, 2, 2)

    And create a custom two dimensional kernel from it:

        >>> gauss_kernel = Model2DKernel(gauss, x_size=9)

    This kernel can now be used like a usual astropy kernel.

    """
    _is_bool = False
    _separable = False

    def __init__(self, model, **kwargs):
        self._separable = False
        if isinstance(model, Fittable2DModel):
            self._model = model
        else:
            raise TypeError("Must be Fittable2DModel")
        super(Model2DKernel, self).__init__(**kwargs)


class PSFKernel(Kernel2D):
    """
    Initialize filter kernel from astropy PSF instance.
    """
    _separable = False

    def __init__(self):
        raise NotImplementedError('Not yet implemented')


class CustomKernel(Kernel):
    """
    Create filter kernel from list or array.

    Parameters
    ----------
    array : list or array
        Filter kernel array. Size must be odd.

    Raises
    ------
    TypeError
        If array is not a list or array.
    KernelSizeError
        If array size is even.

    See also
    --------
    Model2DKernel, Model1DKernel

    Examples
    --------
    Define one dimensional array:

        >>> from astropy.convolution.kernels import CustomKernel
        >>> import numpy as np
        >>> array = np.array([1, 2, 3, 2, 1])
        >>> kernel = CustomKernel(array)
        >>> kernel.dimension
        1

    Define two dimensional array:

        >>> array = np.array([[1, 1, 1], [1, 2, 1], [1, 1, 1]])
        >>> kernel = CustomKernel(array)
        >>> kernel.dimension
        2
    """
    def __init__(self, array):
        self.array = array
        super(CustomKernel, self).__init__(self._array)

    @property
    def array(self):
        """
        Filter kernel array.
        """
        return self._array

    @array.setter
    def array(self, array):
        """
        Filter kernel array setter
        """
        if isinstance(array, np.ndarray):
            self._array = array.astype(np.float64)
        elif isinstance(array, list):
            self._array = np.array(array, dtype=np.float64)
        else:
            raise TypeError("Must be list or array.")

        # Check if array is odd in all axes
        odd = all(axes_size % 2 != 0 for axes_size in self.shape)
        if not odd:
            raise KernelSizeError("Kernel size must be odd in all axes.")

        # Check if array is bool
        ones = self._array == 1.
        zeros = self._array == 0
        self._is_bool = bool(np.all(np.logical_or(ones, zeros)))

        self._truncation = 0.0