summaryrefslogtreecommitdiff
path: root/silx/math/medianfilter/test/test_medianfilter.py
blob: a4c55b2d702f71532c10500a51eb2d6ed9871c63 (plain)
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
# coding: utf-8
# ##########################################################################
# Copyright (C) 2017 European Synchrotron Radiation Facility
#
# Permission is hereby granted, free of charge, to any person obtaining a copy
# of this software and associated documentation files (the "Software"), to deal
# in the Software without restriction, including without limitation the rights
# to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
# copies of the Software, and to permit persons to whom the Software is
# furnished to do so, subject to the following conditions:
#
# The above copyright notice and this permission notice shall be included in
# all copies or substantial portions of the Software.
#
# THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
# IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
# FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
# AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
# LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
# OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
# THE SOFTWARE.
#
# ############################################################################
"""Tests of the median filter"""

__authors__ = ["H. Payno"]
__license__ = "MIT"
__date__ = "02/05/2017"

import unittest
import numpy
from silx.math.medianfilter import medfilt2d
from silx.math.medianfilter.medianfilter import reflect, mirror
from silx.math.medianfilter.medianfilter import MODES as silx_mf_modes
from silx.test.utils import ParametricTestCase
try:
    import scipy
    import scipy.misc
except:
    scipy = None
else:
    import scipy.ndimage

import logging
_logger = logging.getLogger(__name__)

RANDOM_FLOAT_MAT = numpy.array([
    [0.05564293, 0.62717157, 0.75002406, 0.40555336, 0.70278975],
    [0.76532598, 0.02839148, 0.05272484, 0.65166994, 0.42161216],
    [0.23067427, 0.74219128, 0.56049024, 0.44406320, 0.28773158],
    [0.81025249, 0.20303021, 0.68382382, 0.46372299, 0.81281709],
    [0.94691602, 0.07813661, 0.81651256, 0.84220106, 0.33623165]])

RANDOM_INT_MAT = numpy.array([
    [0, 5, 2, 6, 1],
    [2, 3, 1, 7, 1],
    [9, 8, 6, 7, 8],
    [5, 6, 8, 2, 4]])


class TestMedianFilterNearest(ParametricTestCase):
    """Unit tests for the median filter in nearest mode"""

    def testFilter3_100(self):
        """Test median filter on a 10x10 matrix with a 3x3 kernel."""
        dataIn = numpy.arange(100, dtype=numpy.int32)
        dataIn = dataIn.reshape((10, 10))

        dataOut = medfilt2d(image=dataIn,
                            kernel_size=(3, 3),
                            conditional=False,
                            mode='nearest')
        self.assertTrue(dataOut[0, 0] == 1)
        self.assertTrue(dataOut[9, 0] == 90)
        self.assertTrue(dataOut[9, 9] == 98)

        self.assertTrue(dataOut[0, 9] == 9)
        self.assertTrue(dataOut[0, 4] == 5)
        self.assertTrue(dataOut[9, 4] == 93)
        self.assertTrue(dataOut[4, 4] == 44)

    def testFilter3_9(self):
        "Test median filter on a 3x3 matrix with a 3x3 kernel."
        dataIn = numpy.array([0, -1, 1,
                              12, 6, -2,
                              100, 4, 12],
                             dtype=numpy.int16)
        dataIn = dataIn.reshape((3, 3))
        dataOut = medfilt2d(image=dataIn,
                            kernel_size=(3, 3),
                            conditional=False,
                            mode='nearest')
        self.assertTrue(dataOut.shape == dataIn.shape)
        self.assertTrue(dataOut[1, 1] == 4)
        self.assertTrue(dataOut[0, 0] == 0)
        self.assertTrue(dataOut[0, 1] == 0)
        self.assertTrue(dataOut[1, 0] == 6)

    def testFilterWidthOne(self):
        """Make sure a filter of one by one give the same result as the input
        """
        dataIn = numpy.arange(100, dtype=numpy.int32)
        dataIn = dataIn.reshape((10, 10))

        dataOut = medfilt2d(image=dataIn,
                            kernel_size=(1, 1),
                            conditional=False,
                            mode='nearest')

        self.assertTrue(numpy.array_equal(dataIn, dataOut))

    def testFilter3Conditionnal(self):
        """Test that the conditional filter apply correctly in a 10x10 matrix
        with a 3x3 kernel
        """
        dataIn = numpy.arange(100, dtype=numpy.int32)
        dataIn = dataIn.reshape((10, 10))

        dataOut = medfilt2d(image=dataIn,
                            kernel_size=(3, 3),
                            conditional=True,
                            mode='nearest')
        self.assertTrue(dataOut[0, 0] == 1)
        self.assertTrue(dataOut[0, 1] == 1)
        self.assertTrue(numpy.array_equal(dataOut[1:8, 1:8], dataIn[1:8, 1:8]))
        self.assertTrue(dataOut[9, 9] == 98)

    def testFilter3_1D(self):
        """Simple test of a 3x3 median filter on a 1D array"""
        dataIn = numpy.arange(100, dtype=numpy.int32)

        dataOut = medfilt2d(image=dataIn,
                            kernel_size=(5),
                            conditional=False,
                            mode='nearest')

        self.assertTrue(dataOut[0] == 0)
        self.assertTrue(dataOut[9] == 9)
        self.assertTrue(dataOut[99] == 99)


class TestMedianFilterReflect(ParametricTestCase):
    """Unit test for the median filter in reflect mode"""

    def testArange9(self):
        """Test from a 3x3 window to RANDOM_FLOAT_MAT"""
        img = numpy.arange(9, dtype=numpy.int32)
        img = img.reshape(3, 3)
        kernel = (3, 3)
        res = medfilt2d(image=img,
                        kernel_size=kernel,
                        conditional=False,
                        mode='reflect')
        self.assertTrue(
            numpy.array_equal(res.ravel(), [1, 2, 2, 3, 4, 5, 6, 6, 7]))

    def testRandom10(self):
        """Test a (5, 3) window to a RANDOM_FLOAT_MAT"""
        kernel = (5, 3)

        thRes = numpy.array([
            [0.23067427, 0.56049024, 0.56049024, 0.4440632, 0.42161216],
            [0.23067427, 0.62717157, 0.56049024, 0.56049024, 0.46372299],
            [0.62717157, 0.62717157, 0.56049024, 0.56049024, 0.4440632],
            [0.76532598, 0.68382382, 0.56049024, 0.56049024, 0.42161216],
            [0.81025249, 0.68382382, 0.56049024, 0.68382382, 0.46372299]])

        res = medfilt2d(image=RANDOM_FLOAT_MAT,
                        kernel_size=kernel,
                        conditional=False,
                        mode='reflect')

        self.assertTrue(numpy.array_equal(thRes, res))

    def testApplyReflect1D(self):
        """Test the reflect function used for the median filter in reflect mode
        """
        # test for inside values
        self.assertTrue(reflect(2, 3) == 2)
        # test for boundaries values
        self.assertTrue(reflect(3, 3) == 2)
        self.assertTrue(reflect(4, 3) == 1)
        self.assertTrue(reflect(5, 3) == 0)
        self.assertTrue(reflect(6, 3) == 0)
        self.assertTrue(reflect(7, 3) == 1)
        self.assertTrue(reflect(-1, 3) == 0)
        self.assertTrue(reflect(-2, 3) == 1)
        self.assertTrue(reflect(-3, 3) == 2)
        self.assertTrue(reflect(-4, 3) == 2)
        self.assertTrue(reflect(-5, 3) == 1)
        self.assertTrue(reflect(-6, 3) == 0)
        self.assertTrue(reflect(-7, 3) == 0)

    def testRandom10Conditionnal(self):
        """Test the median filter in reflect mode and with the conditionnal
        option"""
        kernel = (3, 1)

        thRes = numpy.array([
            [0.05564293, 0.62717157, 0.75002406, 0.40555336, 0.70278975],
            [0.23067427, 0.62717157, 0.56049024, 0.44406320, 0.42161216],
            [0.76532598, 0.20303021, 0.56049024, 0.46372299, 0.42161216],
            [0.81025249, 0.20303021, 0.68382382, 0.46372299, 0.33623165],
            [0.94691602, 0.07813661, 0.81651256, 0.84220106, 0.33623165]])

        res = medfilt2d(image=RANDOM_FLOAT_MAT,
                        kernel_size=kernel,
                        conditional=True,
                        mode='reflect')
        self.assertTrue(numpy.array_equal(thRes, res))


class TestMedianFilterMirror(ParametricTestCase):
    """Unit test for the median filter in mirror mode
    """

    def testApplyMirror1D(self):
        """Test the reflect function used for the median filter in mirror mode
        """
        # test for inside values
        self.assertTrue(mirror(2, 3) == 2)
        # test for boundaries values
        self.assertTrue(mirror(4, 4) == 2)
        self.assertTrue(mirror(5, 4) == 1)
        self.assertTrue(mirror(6, 4) == 0)
        self.assertTrue(mirror(7, 4) == 1)
        self.assertTrue(mirror(8, 4) == 2)
        self.assertTrue(mirror(-1, 4) == 1)
        self.assertTrue(mirror(-2, 4) == 2)
        self.assertTrue(mirror(-3, 4) == 3)
        self.assertTrue(mirror(-4, 4) == 2)
        self.assertTrue(mirror(-5, 4) == 1)
        self.assertTrue(mirror(-6, 4) == 0)

    def testRandom10(self):
        """Test a (5, 3) window to a random array"""
        kernel = (3, 5)

        thRes = numpy.array([
            [0.05272484, 0.40555336, 0.42161216, 0.42161216, 0.42161216],
            [0.56049024, 0.56049024, 0.4440632, 0.4440632, 0.4440632],
            [0.56049024, 0.46372299, 0.46372299, 0.46372299, 0.46372299],
            [0.68382382, 0.56049024, 0.56049024, 0.46372299, 0.56049024],
            [0.68382382, 0.46372299, 0.68382382, 0.46372299, 0.68382382]])

        res = medfilt2d(image=RANDOM_FLOAT_MAT,
                        kernel_size=kernel,
                        conditional=False,
                        mode='mirror')

        self.assertTrue(numpy.array_equal(thRes, res))

    def testRandom10Conditionnal(self):
        """Test the median filter in reflect mode and with the conditionnal
        option"""
        kernel = (1, 3)

        thRes = numpy.array([
            [0.62717157, 0.62717157, 0.62717157, 0.70278975, 0.40555336],
            [0.02839148, 0.05272484, 0.05272484, 0.42161216, 0.65166994],
            [0.74219128, 0.56049024, 0.56049024, 0.44406320, 0.44406320],
            [0.20303021, 0.68382382, 0.46372299, 0.68382382, 0.46372299],
            [0.07813661, 0.81651256, 0.81651256, 0.81651256, 0.84220106]])

        res = medfilt2d(image=RANDOM_FLOAT_MAT,
                        kernel_size=kernel,
                        conditional=True,
                        mode='mirror')

        self.assertTrue(numpy.array_equal(thRes, res))


class TestMedianFilterShrink(ParametricTestCase):
    """Unit test for the median filter in mirror mode
    """

    def testRandom_3x3(self):
        """Test the median filter in shrink mode and with the conditionnal
        option"""
        kernel = (3, 3)

        thRes = numpy.array([
            [0.62717157, 0.62717157, 0.62717157, 0.65166994, 0.65166994],
            [0.62717157, 0.56049024, 0.56049024, 0.44406320, 0.44406320],
            [0.74219128, 0.56049024, 0.46372299, 0.46372299, 0.46372299],
            [0.74219128, 0.68382382, 0.56049024, 0.56049024, 0.46372299],
            [0.81025249, 0.81025249, 0.68382382, 0.81281709, 0.81281709]])

        res = medfilt2d(image=RANDOM_FLOAT_MAT,
                        kernel_size=kernel,
                        conditional=False,
                        mode='shrink')

        self.assertTrue(numpy.array_equal(thRes, res))

    def testBounds(self):
        """Test the median filter in shrink mode with 3 different kernels
        which should return the same result due to the large values of kernels
        used.
        """
        kernel1 = (1, 9)
        kernel2 = (1, 11)
        kernel3 = (1, 21)

        thRes = numpy.array([[2, 2, 2, 2, 2],
                             [2, 2, 2, 2, 2],
                             [8, 8, 8, 8, 8],
                             [5, 5, 5, 5, 5]])

        resK1 = medfilt2d(image=RANDOM_INT_MAT,
                          kernel_size=kernel1,
                          conditional=False,
                          mode='shrink')

        resK2 = medfilt2d(image=RANDOM_INT_MAT,
                          kernel_size=kernel2,
                          conditional=False,
                          mode='shrink')

        resK3 = medfilt2d(image=RANDOM_INT_MAT,
                          kernel_size=kernel3,
                          conditional=False,
                          mode='shrink')

        self.assertTrue(numpy.array_equal(resK1, thRes))
        self.assertTrue(numpy.array_equal(resK2, resK1))
        self.assertTrue(numpy.array_equal(resK3, resK1))

    def testRandom_3x3Conditionnal(self):
        """Test the median filter in reflect mode and with the conditionnal
        option"""
        kernel = (3, 3)

        thRes = numpy.array([
            [0.05564293, 0.62717157, 0.62717157, 0.40555336, 0.65166994],
            [0.62717157, 0.56049024, 0.05272484, 0.65166994, 0.42161216],
            [0.23067427, 0.74219128, 0.56049024, 0.44406320, 0.46372299],
            [0.81025249, 0.20303021, 0.68382382, 0.46372299, 0.81281709],
            [0.81025249, 0.81025249, 0.81651256, 0.81281709, 0.81281709]])

        res = medfilt2d(image=RANDOM_FLOAT_MAT,
                        kernel_size=kernel,
                        conditional=True,
                        mode='shrink')

        self.assertTrue(numpy.array_equal(res, thRes))

    def testRandomInt(self):
        """Test 3x3 kernel on RANDOM_INT_MAT
        """
        kernel = (3, 3)

        thRes = numpy.array([[3, 2, 5, 2, 6],
                             [5, 3, 6, 6, 7],
                             [6, 6, 6, 6, 7],
                             [8, 8, 7, 7, 7]])

        resK1 = medfilt2d(image=RANDOM_INT_MAT,
                          kernel_size=kernel,
                          conditional=False,
                          mode='shrink')

        self.assertTrue(numpy.array_equal(resK1, thRes))


class TestGeneralExecution(ParametricTestCase):
    """Some general test on median filter application"""

    def testTypes(self):
        """Test that all needed types have their implementation of the median
        filter
        """
        for mode in silx_mf_modes:
            for testType in [numpy.float32, numpy.float64, numpy.int16,
                             numpy.uint16, numpy.int32, numpy.int64,
                             numpy.uint64]:
                with self.subTest(mode=mode, type=testType):
                    data = (numpy.random.rand(10, 10) * 65000).astype(testType)
                    out = medfilt2d(image=data,
                                    kernel_size=(3, 3),
                                    conditional=False,
                                    mode=mode)
                    self.assertTrue(out.dtype.type is testType)

    def testInputDataIsNotModify(self):
        """Make sure input data is not modify by the median filter"""
        dataIn = numpy.arange(100, dtype=numpy.int32)
        dataIn = dataIn.reshape((10, 10))
        dataInCopy = dataIn.copy()

        for mode in silx_mf_modes:
            with self.subTest(mode=mode):
                medfilt2d(image=dataIn,
                          kernel_size=(3, 3),
                          conditional=False,
                          mode=mode)
                self.assertTrue(numpy.array_equal(dataIn, dataInCopy))


def _getScipyAndSilxCommonModes():
    """return the mode which are comparable between silx and scipy"""
    modes = silx_mf_modes.copy()
    del modes['shrink']
    return modes


@unittest.skipUnless(scipy is not None, "scipy not available")
class TestVsScipy(ParametricTestCase):
    """Compare scipy.ndimage.median_filter vs silx.math.medianfilter
    on comparable 
    """
    def testWithArange(self):
        """Test vs scipy with different kernels on arange matrix"""
        data = numpy.arange(10000, dtype=numpy.int32)
        data = data.reshape(100, 100)

        kernels = [(3, 7), (7, 5), (1, 1), (3, 3)]
        modesToTest = _getScipyAndSilxCommonModes()
        for kernel in kernels:
            for mode in modesToTest:
                with self.subTest(kernel=kernel, mode=mode):
                    resScipy = scipy.ndimage.median_filter(input=data,
                                                           size=kernel,
                                                           mode=mode)
                    resSilx = medfilt2d(image=data,
                                        kernel_size=kernel,
                                        conditional=False,
                                        mode=mode)

                    self.assertTrue(numpy.array_equal(resScipy, resSilx))

    def testRandomMatrice(self):
        """Test vs scipy with different kernels on RANDOM_FLOAT_MAT"""
        kernels = [(3, 7), (7, 5), (1, 1), (3, 3)]
        modesToTest = _getScipyAndSilxCommonModes()
        for kernel in kernels:
            for mode in modesToTest:
                with self.subTest(kernel=kernel, mode=mode):
                    resScipy = scipy.ndimage.median_filter(input=RANDOM_FLOAT_MAT,
                                                           size=kernel,
                                                           mode=mode)

                    resSilx = medfilt2d(image=RANDOM_FLOAT_MAT,
                                        kernel_size=kernel,
                                        conditional=False,
                                        mode=mode)

                    self.assertTrue(numpy.array_equal(resScipy, resSilx))

    def testAscentOrLena(self):
        """Test vs scipy with """
        if hasattr(scipy.misc, 'ascent'):
            img = scipy.misc.ascent()
        else:
            img = scipy.misc.lena()

        kernels = [(3, 1), (3, 5), (5, 9), (9, 3)]
        modesToTest = _getScipyAndSilxCommonModes()

        for kernel in kernels:
            for mode in modesToTest:
                with self.subTest(kernel=kernel, mode=mode):
                    resScipy = scipy.ndimage.median_filter(input=img,
                                                           size=kernel,
                                                           mode=mode)

                    resSilx = medfilt2d(image=img,
                                        kernel_size=kernel,
                                        conditional=False,
                                        mode=mode)

                    self.assertTrue(numpy.array_equal(resScipy, resSilx))


def suite():
    test_suite = unittest.TestSuite()
    for test in [TestGeneralExecution, TestVsScipy,
                 TestMedianFilterNearest, TestMedianFilterReflect,
                 TestMedianFilterMirror, TestMedianFilterShrink]:
        test_suite.addTest(
            unittest.defaultTestLoader.loadTestsFromTestCase(test))
    return test_suite


if __name__ == '__main__':
    unittest.main(defaultTest='suite')