-
Notifications
You must be signed in to change notification settings - Fork 49
/
imfilter.jl
1662 lines (1452 loc) · 62.1 KB
/
imfilter.jl
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
# see below for imfilter docstring
# Step 1: if necessary, determine the output's element type
@inline function imfilter(img::AbstractArray, kernel, args...)
imfilter(filter_type(img, kernel), img, kernel, args...)
end
# Step 2: if necessary, put the kernel into cannonical (factored) form
@inline function imfilter(::Type{T}, img::AbstractArray, kernel::Union{ArrayLike,Laplacian}, args...) where T
imfilter(T, img, factorkernel(kernel), args...)
end
@inline function imfilter(::Type{T}, img::AbstractArray{TI}, kernel::AbstractArray{TK}, args...) where {T<:Integer,TI<:Integer,TK<:Integer}
imfilter(T, img, (kernel,), args...)
end
# Step 3: if necessary, fill in the default border
function imfilter(::Type{T}, img::AbstractArray, kernel::ProcessedKernel, args...) where T
imfilter(T, img, kernel, "replicate", args...)
end
function imfilter(::Type{T}, img::AbstractArray, kernel::ProcessedKernel, border::AbstractString, args...) where T
imfilter(T, img, kernel, borderinstance(border), args...)
end
# Step 4: if necessary, allocate the ouput
@inline function imfilter(::Type{T}, img::AbstractArray, kernel::ProcessedKernel, border::AbstractBorder, args...) where T
imfilter!(allocate_output(T, img, kernel, border), img, kernel, border, args...)
end
# Now do the same steps for the case where the user supplies a Resource
@inline function imfilter(r::AbstractResource, img::AbstractArray, kernel, args...)
imfilter(r, filter_type(img, kernel), img, kernel, args...)
end
@inline function imfilter(r::AbstractResource, ::Type{T}, img::AbstractArray, kernel::ArrayLike, args...) where T
imfilter(r, T, img, factorkernel(kernel), args...)
end
# For steps 3 & 4, we make args... explicit as a means to prevent
# specifying both r and an algorithm
function imfilter(r::AbstractResource, ::Type{T}, img::AbstractArray, kernel::ProcessedKernel) where T
imfilter(r, T, img, kernel, Pad(:replicate)) # supply the default border
end
function imfilter(r::AbstractResource, ::Type{T}, img::AbstractArray, kernel::ProcessedKernel, border::AbstractString) where T
imfilter(r, T, img, kernel, borderinstance(border))
end
function imfilter(r::AbstractResource, ::Type{T}, img::AbstractArray, kernel::ProcessedKernel, border::AbstractBorder) where T
imfilter!(r, allocate_output(T, img, kernel, border), img, kernel, border)
end
"""
imfilter([T], img, kernel, [border="replicate"], [alg]) --> imgfilt
imfilter([r], img, kernel, [border="replicate"], [alg]) --> imgfilt
imfilter(r, T, img, kernel, [border="replicate"], [alg]) --> imgfilt
Filter a one, two or multidimensional array `img` with a `kernel` by computing
their correlation.
# Details
The term *filtering* emerges in the context of a Fourier transformation of
an image, which maps an image from its canonical spatial domain to its
concomitant frequency domain. Manipulating an image in the frequency domain
amounts to retaining or discarding particular frequency components—a process
analogous to sifting or filtering [1]. Because the Fourier transform establishes a
link between the spatial and frequency representation of an image, one can
interpret various image manipulations in the spatial domain as filtering
operations which accept or reject specific frequencies.
The phrase *spatial filtering* is often used to emphasise that an operation
is, at least conceptually, devised in the context of the spatial domain of an
image. One further distinguishes between linear and non-linear spatial
filtering. A filter is called linear if the operation performed on the pixels is
linear, and is labeled non-linear otherwise.
An image filter can be represented by a function
```math
w: \\{s\\in \\mathbb{Z} \\mid -k_1 \\le s \\le k_1 \\} \\times \\{t \\in \\mathbb{Z} \\mid -k_2 \\le t \\le k_2 \\} \\rightarrow \\mathbb{R},
```
where ``k_i \\in \\mathbb{N}`` (i = 1,2). It is common to define ``k_1 = 2a+1``
and ``k_2 = 2b + 1``, where ``a`` and ``b`` are integers, which ensures that the
filter dimensions are of odd size. Typically, ``k_1`` equals ``k_2`` and so,
dropping the subscripts, one speaks of a ``k \\times k`` filter. Since the
domain of the filter represents a grid of spatial coordinates, the filter is
often called a mask and is visualized as a grid. For example, a ``3 \\times 3``
mask can be potrayed as follows:
```math
\\scriptsize
\\begin{matrix}
\\boxed{
\\begin{matrix}
\\phantom{w(-9,-9)} \\\\
w(-1,-1) \\\\
\\phantom{w(-9,-9)} \\\\
\\end{matrix}
}
&
\\boxed{
\\begin{matrix}
\\phantom{w(-9,-9)} \\\\
w(-1,0) \\\\
\\phantom{w(-9,-9)} \\\\
\\end{matrix}
}
&
\\boxed{
\\begin{matrix}
\\phantom{w(-9,-9)} \\\\
w(-1,1) \\\\
\\phantom{w(-9,-9)} \\\\
\\end{matrix}
}
\\\\
\\\\
\\boxed{
\\begin{matrix}
\\phantom{w(-9,-9)} \\\\
w(0,-1) \\\\
\\phantom{w(-9,-9)} \\\\
\\end{matrix}
}
&
\\boxed{
\\begin{matrix}
\\phantom{w(-9,-9)} \\\\
w(0,0) \\\\
\\phantom{w(-9,-9)} \\\\
\\end{matrix}
}
&
\\boxed{
\\begin{matrix}
\\phantom{w(-9,-9)} \\\\
w(0,1) \\\\
\\phantom{w(-9,-9)} \\\\
\\end{matrix}
}
\\\\
\\\\
\\boxed{
\\begin{matrix}
\\phantom{w(-9,-9)} \\\\
w(1,-1) \\\\
\\phantom{w(-9,-9)} \\\\
\\end{matrix}
}
&
\\boxed{
\\begin{matrix}
\\phantom{w(-9,-9)} \\\\
w(1,0) \\\\
\\phantom{w(-9,-9)} \\\\
\\end{matrix}
}
&
\\boxed{
\\begin{matrix}
\\phantom{w(-9,-9)} \\\\
w(1,1) \\\\
\\phantom{w(-9,-9)} \\\\
\\end{matrix}
}
\\end{matrix}.
```
The values of ``w(s,t)`` are referred to as *filter coefficients*.
## Discrete convolution versus correlation
There are two fundamental and closely related operations that one regularly
performs on an image with a filter. The operations are called discrete
*correlation* and *convolution*.
The correlation operation, denoted by the symbol ``\\star``, is given in two
dimensions by the expression
```math
\\begin{aligned}
g(x,y) = w(x,y) \\star f(x,y) = \\sum_{s = -a}^{a} \\sum_{t=-b}^{b} w(s,t) f(x+s, y+t),
\\end{aligned}
```
whereas the comparable convolution operation, denoted by the symbol ``\\ast``,
is given in two dimensions by
```math
\\begin{aligned}
h(x,y) = w(x,y) \\ast f(x,y) = \\sum_{s = -a}^{a} \\sum_{t=-b}^{b} w(s,t) f(x-s, y-t).
\\end{aligned}
```
Since a digital image is of finite extent, both of these operations are
undefined at the borders of the image. In particular, for an image of size ``M
\\times N``, the function ``f(x \\pm s, y \\pm t)`` is only defined for ``1 \\le
x \\pm s \\le N`` and ``1 \\le y \\pm t \\le M``. In practice one addresses this
problem by artificially expanding the domain of the image. For example, one can
pad the image with zeros. Other padding strategies are possible, and they are
discussed in more detail in the *Options* section of this documentation.
## One-dimensional illustration
The difference between correlation and convolution is best understood with
recourse to a one-dimensional example adapted from [1]. Suppose that a filter
``w:\\{-1,0,1\\}\\rightarrow \\mathbb{R}`` has coefficients
```math
\\begin{matrix}
\\boxed{1} & \\boxed{2} & \\boxed{3}
\\end{matrix}.
```
Consider a discrete unit impulse function ``f: \\{x \\in \\mathbb{Z} \\mid 1 \\le x
\\le 7 \\} \\rightarrow \\{0,1\\}`` that has been padded with zeros. The function
can be visualised as an image
```math
\\boxed{
\\begin{matrix}
0 & \\boxed{0} & \\boxed{0} & \\boxed{0} & \\boxed{1} & \\boxed{0} & \\boxed{0} & \\boxed{0} & 0
\\end{matrix}}.
```
The correlation operation can be interpreted as sliding ``w`` along the image
and computing the sum of products at each location. For example,
```math
\\begin{matrix}
0 & 0 & 0 & 0 & 1 & 0 & 0 & 0 & 0 \\\\
1 & 2 & 3 & & & & & & \\\\
& 1 & 2 & 3 & & & & & \\\\
& & 1 & 2 & 3 & & & & \\\\
& & & 1 & 2 & 3 & & & \\\\
& & & & 1 & 2 & 3 & & \\\\
& & & & & 1 & 2 & 3 & \\\\
& & & & & & 1 & 2 & 3,
\\end{matrix}
```
yields the output ``g: \\{x \\in \\mathbb{Z} \\mid 1 \\le x \\le 7 \\} \\rightarrow
\\mathbb{R}``, which when visualized as a digital image, is equal to
```math
\\boxed{
\\begin{matrix}
\\boxed{0} & \\boxed{0} & \\boxed{3} & \\boxed{2} & \\boxed{1} & \\boxed{0} & \\boxed{0}
\\end{matrix}}.
```
The interpretation of the convolution operation is analogous to correlation,
except that the filter ``w`` has been rotated by 180 degrees. In particular,
```math
\\begin{matrix}
0 & 0 & 0 & 0 & 1 & 0 & 0 & 0 & 0 \\\\
3 & 2 & 1 & & & & & & \\\\
& 3 & 2 & 1 & & & & & \\\\
& & 3 & 2 & 1 & & & & \\\\
& & & 3 & 2 & 1 & & & \\\\
& & & & 3 & 2 & 1 & & \\\\
& & & & & 3 & 2 & 1 & \\\\
& & & & & & 3 & 2 & 1,
\\end{matrix}
```
yields the output ``h: \\{x \\in \\mathbb{Z} \\mid 1 \\le x \\le 7 \\} \\rightarrow \\mathbb{R}`` equal to
```math
\\boxed{
\\begin{matrix}
\\boxed{0} & \\boxed{0} & \\boxed{1} & \\boxed{2} & \\boxed{3} & \\boxed{0} & \\boxed{0}
\\end{matrix}}.
```
Instead of rotating the filter mask, one could instead rotate ``f`` and still
obtained the same convolution result. In fact, the conventional notation for
convolution indicates that ``f`` is flipped and not ``w``. If ``w`` is
symmetric, then convolution and correlation give the same outcome.
### Two-dimensional illustration
For a two-dimensional example, suppose the filter ``w:\\{-1, 0 ,1\\} \\times
\\{-1,0,1\\} \\rightarrow \\mathbb{R}`` has coefficients
```math
\\begin{matrix}
\\boxed{1} & \\boxed{2} & \\boxed{3} \\\\ \\\\
\\boxed{4} & \\boxed{5} & \\boxed{6} \\\\ \\\\
\\boxed{7} & \\boxed{8} & \\boxed{9}
\\end{matrix},
```
and consider a two-dimensional discrete unit impulse function
```math
f:\\{x \\in \\mathbb{Z} \\mid 1 \\le x \\le 7 \\} \\times \\{y \\in \\mathbb{Z} \\mid 1 \\le y \\le 7 \\}\\rightarrow \\{ 0,1\\}
```
that has been padded with zeros:
```math
\\boxed{
\\begin{matrix}
0 & 0 & 0 & 0 & 0 & 0 & 0 \\\\ \\\\
0 & \\boxed{0} & \\boxed{0} & \\boxed{0} & \\boxed{0} & \\boxed{0} & 0 \\\\ \\\\
0 & \\boxed{0} & \\boxed{0} & \\boxed{0} & \\boxed{0} & \\boxed{0} & 0 \\\\ \\\\
0 & \\boxed{0} & \\boxed{0} & \\boxed{1} & \\boxed{0} & \\boxed{0} & 0 \\\\ \\\\
0 & \\boxed{0} & \\boxed{0} & \\boxed{0} & \\boxed{0} & \\boxed{0} & 0 \\\\ \\\\
0 & \\boxed{0} & \\boxed{0} & \\boxed{0} & \\boxed{0} & \\boxed{0} & 0 \\\\ \\\\
0 & 0 & 0 & 0 & 0 & 0 & 0
\\end{matrix}}.
```
The correlation operation ``w(x,y) \\star f(x,y)`` yields the output
```math
\\boxed{
\\begin{matrix}
\\boxed{0} & \\boxed{0} & \\boxed{0} & \\boxed{0} & \\boxed{0} \\\\ \\\\
\\boxed{0} & \\boxed{9} & \\boxed{8} & \\boxed{7} & \\boxed{0} \\\\ \\\\
\\boxed{0} & \\boxed{6} & \\boxed{5} & \\boxed{4} & \\boxed{0} \\\\ \\\\
\\boxed{0} & \\boxed{3} & \\boxed{2} & \\boxed{1} & \\boxed{0} \\\\ \\\\
\\boxed{0} & \\boxed{0} & \\boxed{0} & \\boxed{0} & \\boxed{0}
\\end{matrix}},
```
whereas the convolution operation ``w(x,y) \\ast f(x,y)`` produces
```math
\\boxed{
\\begin{matrix}
\\boxed{0} & \\boxed{0} & \\boxed{0} & \\boxed{0} & \\boxed{0} \\\\ \\\\
\\boxed{0} & \\boxed{1} & \\boxed{2} & \\boxed{3} & \\boxed{0}\\\\ \\\\
\\boxed{0} & \\boxed{4} & \\boxed{5} & \\boxed{6} & \\boxed{0} \\\\ \\\\
\\boxed{0} & \\boxed{7} & \\boxed{8} & \\boxed{9} & \\boxed{0} \\\\ \\\\
\\boxed{0} & \\boxed{0} & \\boxed{0} & \\boxed{0} & \\boxed{0}
\\end{matrix}}.
```
## Discrete convolution and correlation as matrix multiplication
Discrete convolution and correlation operations can also be formulated as a
matrix multiplication, where one of the inputs is converted to a [Toeplitz](https://en.wikipedia.org/wiki/Toeplitz_matrix)
matrix, and the other is represented as a column vector. For example, consider a
function ``f:\\{x \\in \\mathbb{N} \\mid 1 \\le x \\le M \\} \\rightarrow \\mathbb{R}``
and a filter ``w: \\{s \\in \\mathbb{N} \\mid -k_1 \\le s \\le k_1 \\} \\rightarrow
\\mathbb{R}``. Then the matrix multiplication
```math
\\begin{bmatrix}
w(-k_1) & 0 & \\ldots & 0 & 0 \\\\
\\vdots & w(-k_1) & \\ldots & \\vdots & 0 \\\\
w(k_1) & \\vdots & \\ldots & 0 & \\vdots \\\\
0 & w(k_1) & \\ldots & w(-k_1) & 0 \\\\
0 & 0 & \\ldots & \\vdots & w(-k_1) \\\\
\\vdots & \\vdots & \\ldots & w(k_1) & \\vdots \\\\
0 & 0 & 0 & 0 & w(k_1)
\\end{bmatrix}
\\begin{bmatrix}
f(1) \\\\
f(2) \\\\
f(3) \\\\
\\vdots \\\\
f(M)
\\end{bmatrix}
```
is equivalent to the convolution ``w(s) \\ast f(x)`` assuming that the border of
``f(x)`` has been padded with zeros.
To represent multidimensional convolution as matrix multiplication one
reshapes the multidimensional arrays into column vectors and proceeds in an
analogous manner. Naturally, the result of the matrix multiplication will need
to be reshaped into an appropriate multidimensional array.
# Options
The following subsections describe valid options for the function arguments in
more detail.
## Choices for `r`
You can dispatch to different implementations by passing in a resource `r`
as defined by the [ComputationalResources](https://github.com/timholy/ComputationalResources.jl) package.
For example,
```julia
imfilter(ArrayFireLibs(), img, kernel)
```
would request that the computation be performed on the GPU using the
ArrayFire libraries.
## Choices for `T`
Optionally, you can control the element type of the output image by
passing in a type `T` as the first argument.
## Choices for `img`
You can specify a one, two or multidimensional array defining your image.
## Choices for `kernel`
The `kernel[0,0,..]` parameter corresponds to the origin (zero displacement) of
the kernel; you can use `centered` to place the origin at the array center, or
use the OffsetArrays package to set `kernel`'s indices manually. For example, to
filter with a random *centered* 3x3 kernel, you could use either of the
following:
kernel = centered(rand(3,3))
kernel = OffsetArray(rand(3,3), -1:1, -1:1)
The `kernel` parameter can be specified as an array or as a "factored kernel", a
tuple `(filt1, filt2, ...)` of filters to apply along each axis of the image. In
cases where you know your kernel is separable, this format can speed processing.
Each of these should have the same dimensionality as the image itself, and be
shaped in a manner that indicates the filtering axis, e.g., a 3x1 filter for
filtering the first dimension and a 1x3 filter for filtering the second
dimension. In two dimensions, any kernel passed as a single matrix is checked
for separability; if you want to eliminate that check, pass the kernel as a
single-element tuple, `(kernel,)`.
## Choices for `border`
At the image edge, `border` is used to specify the padding which will be used
to extrapolate the image beyond its original bounds. As an indicative example
of each option the results of the padding are illustrated on an image consisting of
a row of six pixels which are specified alphabetically: ``\\boxed{a \\, b \\, c \\, d \\, e \\, f}``.
We show the effects of padding only on the left and right border, but analogous
consequences hold for the top and bottom border.
### `"replicate"` (default)
The border pixels extend beyond the image boundaries.
```math
\\boxed{
\\begin{array}{l|c|r}
a\\, a\\, a\\, a & a \\, b \\, c \\, d \\, e \\, f & f \\, f \\, f \\, f
\\end{array}
}
```
See also: [`Pad`](@ref), [`padarray`](@ref), [`Inner`](@ref), [`NA`](@ref) and
[`NoPad`](@ref)
### `"circular"`
The border pixels wrap around. For instance, indexing beyond the left border
returns values starting from the right border.
```math
\\boxed{
\\begin{array}{l|c|r}
c\\, d\\, e\\, f & a \\, b \\, c \\, d \\, e \\, f & a \\, b \\, c \\, d
\\end{array}
}
```
See also: [`Pad`](@ref), [`padarray`](@ref), [`Inner`](@ref), [`NA`](@ref) and
[`NoPad`](@ref)
### `"reflect"`
The border pixels reflect relative to a position between pixels. That is, the
border pixel is omitted when mirroring.
```math
\\boxed{
\\begin{array}{l|c|r}
e\\, d\\, c\\, b & a \\, b \\, c \\, d \\, e \\, f & e \\, d \\, c \\, b
\\end{array}
}
```
See also: [`Pad`](@ref), [`padarray`](@ref), [`Inner`](@ref), [`NA`](@ref) and
[`NoPad`](@ref)
### `"symmetric"`
The border pixels reflect relative to the edge itself.
```math
\\boxed{
\\begin{array}{l|c|r}
d\\, c\\, b\\, a & a \\, b \\, c \\, d \\, e \\, f & f \\, e \\, d \\, c
\\end{array}
}
```
See also: [`Pad`](@ref), [`padarray`](@ref), [`Inner`](@ref), [`NA`](@ref) and
[`NoPad`](@ref)
### `Fill(m)`
The border pixels are filled with a specified value ``m``.
```math
\\boxed{
\\begin{array}{l|c|r}
m\\, m\\, m\\, m & a \\, b \\, c \\, d \\, e \\, f & m \\, m \\, m \\, m
\\end{array}
}
```
See also: [`Pad`](@ref), [`padarray`](@ref), [`Inner`](@ref), [`NA`](@ref) and
[`NoPad`](@ref)
### `Inner()`
Indicate that edges are to be discarded in filtering, only the interior of the
result is to be returned.
See also: [`Pad`](@ref), [`padarray`](@ref), [`Inner`](@ref), [`NA`](@ref) and
[`NoPad`](@ref)
### `NA()`
Choose filtering using "NA" (Not Available) boundary conditions. This
is most appropriate for filters that have only positive weights, such
as blurring filters.
See also: [`Pad`](@ref), [`padarray`](@ref), [`Inner`](@ref), [`NA`](@ref) and
[`NoPad`](@ref)
## Choices for `alg`
The `alg` parameter allows you to choose the particular algorithm: `FIR()`
(finite impulse response, aka traditional digital filtering) or `FFT()`
(Fourier-based filtering). If no choice is specified, one will be chosen based
on the size of the image and kernel in a way that strives to deliver good
performance. Alternatively you can use a custom filter type, like
[`KernelFactors.IIRGaussian`](@ref).
# Examples
The following subsections highlight some common use cases.
## Convolution versus correlation
```julia
# Create a two-dimensional discrete unit impulse function.
f = fill(0,(9,9));
f[5,5] = 1;
# Specify a filter coefficient mask and set the center of the mask as the origin.
w = centered([1 2 3; 4 5 6 ; 7 8 9]);
#=
The default operation of `imfilter` is correlation. By reflecting `w` we
compute the convolution of `f` and `w`. `Fill(0,w)` indicates that we wish to
pad the border of `f` with zeros. The amount of padding is automatically
determined by considering the length of w.
=#
correlation = imfilter(f,w,Fill(0,w))
convolution = imfilter(f,reflect(w),Fill(0,w))
```
## Miscellaneous border padding options
```julia
# Example function values f, and filter coefficients w.
f = reshape(1.0:81.0,9,9)
w = centered(reshape(1.0:9.0,3,3))
# You can designate the type of padding by specifying an appropriate string.
imfilter(f,w,"replicate")
imfilter(f,w,"circular")
imfilter(f,w,"symmetric")
imfilter(f,w,"reflect")
# Alternatively, you can explicitly use the Pad type to designate the padding style.
imfilter(f,w,Pad(:replicate))
imfilter(f,w,Pad(:circular))
imfilter(f,w,Pad(:symmetric))
imfilter(f,w,Pad(:reflect))
# If you want to pad with a specific value then use the Fill type.
imfilter(f,w,Fill(0,w))
imfilter(f,w,Fill(1,w))
imfilter(f,w,Fill(-1,w))
#=
Specify 'Inner()' if you want to retrieve the interior sub-array of f for which
the filtering operation is defined without padding.
=#
imfilter(f,w,Inner())
```
# References
1. R. C. Gonzalez and R. E. Woods. *Digital Image Processing (3rd Edition)*. Upper Saddle River, NJ, USA: Prentice-Hall, 2006.
See also: [`imfilter!`](@ref), [`centered`](@ref), [`padarray`](@ref), [`Pad`](@ref), [`Fill`](@ref), [`Inner`](@ref), [`KernelFactors.IIRGaussian`](@ref).
"""
imfilter
# see below for imfilter! docstring
# imfilter! can be called directly, so we take steps 2&3 here too. We
# have to be a little more cautious to make sure that later methods
# don't inadvertently call back to these: in methods that take an
# AbstractResource argument, exclude `NoPad()` as a border option.
function imfilter!(out::AbstractArray, img::AbstractArray, kernel::Union{ArrayLike,Laplacian}, args...)
imfilter!(out, img, factorkernel(kernel), args...)
end
function imfilter!(r::AbstractResource, out::AbstractArray, img::AbstractArray, kernel, args...)
imfilter!(r, out, img, factorkernel(kernel), args...)
end
function imfilter!(out::AbstractArray, img::AbstractArray, kernel::ProcessedKernel)
imfilter!(out, img, kernel, Pad(:replicate))
end
function imfilter!(r::AbstractResource, out::AbstractArray, img::AbstractArray, kernel::ProcessedKernel)
imfilter!(r, out, img, kernel, Pad(:replicate))
end
function imfilter!(out::AbstractArray, img::AbstractArray, kernel::ProcessedKernel, border::AbstractString, args...)
imfilter!(out, img, kernel, borderinstance(border), args...)
end
function imfilter!(r::AbstractResource, out::AbstractArray, img::AbstractArray, kernel::ProcessedKernel, border::AbstractString)
imfilter!(r, out, img, kernel, borderinstance(border))
end
# Step 5: if necessary, pick an algorithm
function imfilter!(out::AbstractArray, img::AbstractArray, kernel::ProcessedKernel, border::AbstractBorder)
imfilter!(out, img, kernel, border, filter_algorithm(out, img, kernel))
end
function imfilter!(out::AbstractArray, img::AbstractArray, kernel::ProcessedKernel, border::AbstractBorder, alg::Alg)
local ret
try
ret = imfilter!(default_resource(alg_defaults(alg, out, kernel)), out, img, kernel, border)
catch err
if isa(err, InexactError)
Tw = Float64
if eltype(img) <: Integer
try
# If a type doesn't support widen, it would be bad
# if our attempt to be helpful triggered a
# completely different error...
Tw = widen(eltype(img))
catch
end
end
@warn "Likely overflow or conversion error detected. Consider specifying the output type, e.g., `imfilter($Tw, img, kernel, ...)`"
end
rethrow(err)
end
ret
end
"""
imfilter!(imgfilt, img, kernel, [border="replicate"], [alg])
imfilter!(r, imgfilt, img, kernel, border, [inds])
imfilter!(r, imgfilt, img, kernel, border::NoPad, [inds=axes(imgfilt)])
Filter an array `img` with kernel `kernel` by computing their
correlation, storing the result in `imgfilt`.
The indices of `imgfilt` determine the region over which the filtered
image is computed---you can use this fact to select just a specific
region of interest, although be aware that the input `img` might still
get padded. Alteratively, explicitly provide the indices `inds` of
`imgfilt` that you want to calculate, and use `NoPad` boundary
conditions. In such cases, you are responsible for supplying
appropriate padding: `img` must be indexable for all of the locations
needed for calculating the output. This syntax is best-supported for
FIR filtering; in particular, that that IIR filtering can lead to
results that are inconsistent with respect to filtering the entire
array.
See also: [`imfilter`](@ref).
"""
imfilter!
# Step 6: pad the input
# NA "padding": normalizing by the number of available values (similar to nanmean)
function imfilter!(r::AbstractResource,
out::AbstractArray{S,N},
img::AbstractArray{T,N},
kernel::ProcessedKernel,
border::NA{0}) where {T,S,N}
_imfilter_na!(r, out, img, kernel, border)
end
function _imfilter_na!(r::AbstractResource,
out::AbstractArray{S,N},
img::AbstractArray{T,N},
kernel::ProcessedKernel,
border::NA{0}) where {T,S,N}
nanflag = isnan.(img)
hasnans = any(nanflag)
if hasnans || !isseparable(kernel)
imfilter_na_inseparable!(r, out, img, nanflag, kernel)
else
imfilter_na_separable!(r, out, img, kernel)
end
out
end
# for types that can't have NaNs, we can skip the isnan check
function _imfilter_na!(r::AbstractResource,
out::AbstractArray{S,N},
img::AbstractArray{T,N},
kernel::ProcessedKernel,
border::NA{0}) where {S,T<:Union{Integer,FixedColorant},N}
if isseparable(kernel)
imfilter_na_separable!(r, out, img, kernel)
else
nanflag = fill(false, axes(img))
imfilter_na_inseparable!(r, out, img, nanflag, kernel)
end
end
# Any other kind of not-fully-specified padding
function imfilter!(r::AbstractResource,
out::AbstractArray{S,N},
img::AbstractArray{T,N},
kernel::ProcessedKernel,
border::BorderSpecAny) where {S,T,N}
bord = border(kernel, img, Alg(r)) # if it's FFT, the size of img is also relevant
imfilter!(r, out, img, kernel, bord)
end
# Any fully-specified padding
function imfilter!(r::AbstractResource,
out::AbstractArray{S,N},
img::AbstractArray{T,N},
kernel::ProcessedKernel,
border::AbstractBorder) where {S,T,N}
A = padarray(S, img, border)
# By specifying NoPad(), we ensure that dispatch will never
# accidentally "go back" to an earlier routine and apply more
# padding
imfilter!(r, out, A, kernel, NoPad(border))
end
# # An optimized case that performs only "virtual padding"
# function imfilter!{S,T,N,A<:Union{FIR,FIRTiled}}(r::AbstractCPU{A},
# out::AbstractArray{S,N},
# img::AbstractArray{T,N},
# kernel::ProcessedKernel,
# border::Pad{0})
# # The fast path: handle the points that don't need padding
# iinds = map(intersect, interior(img, kernel), axes(out))
# imfilter!(r, out, img, kernel, NoPad(border), iinds)
# # The not-so-fast path: handle the edges
# # TODO: when the kernel is factored, move this logic in to each factor
# # This is especially important for bigger kernels, where the product pkernel is larger
# padded = view(img, padindices(img, border(kernel))...)
# pkernel = kernelconv(kernel...)
# _imfilter_iter!(r, out, padded, pkernel, EdgeIterator(axes(out), iinds))
# end
### "Scheduler" methods (all with NoPad)
# These methods handle much of what Halide calls "the schedule."
# Together they handle the order-of-operations for separable and/or
# cascaded kernels, and even implement multithreadable tiling for FIR
# filtering.
# Trivial kernel (a copy operation)
function imfilter!(r::AbstractResource, out::AbstractArray, A::AbstractArray, kernel::Tuple{}, ::NoPad, inds::Indices=axes(out))
R = CartesianIndices(inds)
copyto!(out, R, A, R)
end
# A single kernel
function imfilter!(r::AbstractResource, out::AbstractArray, A::AbstractArray, kernel::Tuple{Any}, border::NoPad, inds::Indices=axes(out))
kern = kernel[1]
iscopy(kern) && return imfilter!(r, out, A, (), border, inds)
imfilter!(r, out, A, samedims(out, kern), border, inds)
end
# A filter cascade (2 or more filters)
function imfilter!(r::AbstractResource, out::AbstractArray, A::AbstractArray, kernel::Tuple{Any,Any,Vararg{Any}}, border::NoPad, inds=axes(out))
kern = kernel[1]
iscopy(kern) && return imfilter!(r, out, A, tail(kernel), border, inds)
# For multiple stages of filtering, we introduce a second buffer
# and swap them at each stage. The first of the two is the one
# that holds the most recent result.
A2 = tempbuffer(A, eltype(out), kernel)
indsstep = shrink(expand(inds, calculate_padding(kernel)), kern)
_imfilter!(r, out, A, A2, kernel, border, indsstep)
return out
end
### Use a tiled algorithm for the cascaded case
function imfilter!(r::AbstractCPU{FIRTiled{N}}, out::AbstractArray{S,N}, A::AbstractArray{T,N}, kernel::Tuple{Any,Any,Vararg{Any}}, border::NoPad, inds=axes(out)) where {S,T,N}
kern = kernel[1]
iscopy(kern) && return imfilter!(r, out, A, tail(kernel), border, inds)
tmp = tile_allocate(filter_type(A, kernel), r.settings.tilesize, kernel)
_imfilter_tiled!(r, out, A, kernel, border, tmp, inds)
out
end
### Scheduler support methods
## No tiling, filter cascade
# For these (internal) methods, `indsstep` refers to `inds` for this
# step of filtering, not the indices of `out` that we want to finally
# target.
# When `kernel` is (originally) a tuple that has both TriggsSdika and
# FIR filters, the overall padding gets doubled, yet we only trim off
# the minimum at each stage. Consequently, `indsstep` might be
# optimistic about the range available in `out`; therefore we use
# `intersect`.
function _imfilter!(r, out::AbstractArray, A1, A2, kernel::Tuple{}, border::NoPad, indsstep::Indices)
imfilter!(r, out, A1, kernel, border, map(intersect, indsstep, axes(out)))
end
function _imfilter!(r, out::AbstractArray, A1, A2, kernel::Tuple{Any}, border::NoPad, indsstep::Indices)
imfilter!(r, out, A1, kernel[1], border, map(intersect, indsstep, axes(out)))
end
# For IIR, it's important to filter over the whole passed-in range,
# and then copy! to out
function _imfilter!(r, out::AbstractArray, A1, A2, kernel::Tuple{AnyIIR}, border::NoPad, indsstep::Indices)
if indsstep != axes(out)
imfilter!(r, A2, A1, kernel[1], border, indsstep)
R = CartesianIndices(map(intersect, indsstep, axes(out)))
return copyto!(out, R, A2, R)
end
imfilter!(r, out, A1, kernel[1], border, indsstep)
end
function _imfilter!(r, out::AbstractArray, A1, A2::AbstractArray, kernel::Tuple{Any,Any,Vararg{Any}}, border::NoPad, indsstep::Indices)
kern = kernel[1]
iscopy(kern) && return _imfilter!(r, out, A1, A2, tail(kernel), border, indsstep)
kernN = samedims(A2, kern)
imfilter!(r, A2, A1, kernN, border, indsstep) # store result in A2
kernelt = tail(kernel)
newinds = next_shrink(indsstep, kernelt)
_imfilter!(r, out, A2, A1, tail(kernel), border, newinds) # swap the buffers
end
function _imfilter!(r, out::AbstractArray, A1, A2::Tuple{AbstractArray,AbstractArray}, kernel::Tuple{Any,Any,Vararg{Any}}, border::NoPad, indsstep::Indices)
kern = kernel[1]
iscopy(kern) && return _imfilter!(r, out, A1, A2, tail(kernel), border, indsstep)
A2_1, A2_2 = A2
kernN = samedims(A2_1, kern)
imfilter!(r, A2_1, A1, kernN, border, indsstep) # store result in A2
kernelt = tail(kernel)
newinds = next_shrink(indsstep, kernelt)
_imfilter!(r, out, A2_1, A2_2, tail(kernel), border, newinds) # swap the buffers
end
# Single-threaded, pair of kernels (with only one temporary buffer required)
function _imfilter_tiled!(r::CPU1, out, A, kernel::Tuple{Any,Any}, border::NoPad, tiles::Vector{AA}, indsout) where AA<:AbstractArray
k1, k2 = kernel
tile = tiles[1]
indsk2, indstile = axes(k2), axes(tile)
sz = map(length, indstile)
chunksz = map(length, shrink(indstile, indsk2))
for tinds in TileIterator(indsout, chunksz)
tileinds = expand(tinds, k2)
tileb = TileBuffer(tile, tileinds)
imfilter!(r, tileb.view, A, samedims(tileb, k1), border, tileinds)
imfilter!(r, out, tileb.view, samedims(out, k2), border, tinds)
end
out
end
# Multithreaded, pair of kernels
function _imfilter_tiled!(r::CPUThreads, out, A, kernel::Tuple{Any,Any}, border::NoPad, tiles::Vector{AA}, indsout) where AA<:AbstractArray
k1, k2 = kernel
tile = tiles[1]
indsk2, indstile = axes(k2), axes(tile)
sz = map(length, indstile)
chunksz = map(length, shrink(indstile, indsk2))
tileinds_all = collect(expand(inds, k2) for inds in TileIterator(indsout, chunksz))
_imfilter_tiled_threads!(CPU1(r), out, A, samedims(out, k1), samedims(out, k2), border, tileinds_all, tiles)
end
# This must be in a separate function due to #15276
@noinline function _imfilter_tiled_threads!(r1, out, A, k1, k2, border, tileinds_all, tile::Vector{AA}) where AA<:AbstractArray
Threads.@threads for i = 1:length(tileinds_all)
id = Threads.threadid()
tileinds = tileinds_all[i]
tileb = TileBuffer(tile[id], tileinds)
imfilter!(r1, tileb, A, k1, border, tileinds)
imfilter!(r1, out, tileb, k2, border, shrink(tileinds, k2))
end
out
end
# Single-threaded, multiple kernels (requires two tile buffers, swapping on each iteration)
function _imfilter_tiled!(r::CPU1, out, A, kernel::Tuple{Any,Any,Vararg{Any}}, border::NoPad, tiles::Vector{Tuple{AA,AA}}, indsout) where AA<:AbstractArray
k1, kt = kernel[1], tail(kernel)
tilepair = tiles[1]
indstile = axes(tilepair[1])
sz = map(length, indstile)
chunksz = map(length, shrink(indstile, kt))
for tinds in TileIterator(indsout, chunksz)
tileinds = expand(tinds, kt)
tileb1 = TileBuffer(tilepair[1], tileinds)
imfilter!(r, tileb1, A, samedims(tileb1, k1), border, tileinds)
_imfilter_tiled_swap!(r, out, kt, border, (tileb1, tilepair[2]))
end
end
# Multithreaded, multiple kernels
function _imfilter_tiled!(r::CPUThreads, out, A, kernel::Tuple{Any,Any,Vararg{Any}}, border::NoPad, tiles::Vector{Tuple{AA,AA}}, indsout) where AA<:AbstractArray
k1, kt = kernel[1], tail(kernel)
tilepair = tiles[1]
indstile = axes(tilepair[1])
sz = map(length, indstile)
chunksz = map(length, shrink(indstile, kt))
tileinds_all = collect(expand(inds, kt) for inds in TileIterator(indsout, chunksz))
_imfilter_tiled_threads!(CPU1(r), out, A, samedims(out, k1), kt, border, tileinds_all, tiles)
end
# This must be in a separate function due to #15276
@noinline function _imfilter_tiled_threads!(r1, out, A, k1, kt, border, tileinds_all, tiles::Vector{Tuple{AA,AA}}) where AA<:AbstractArray
Threads.@threads for i = 1:length(tileinds_all)
tileinds = tileinds_all[i]
id = Threads.threadid()
tile1, tile2 = tiles[id]
tileb1 = TileBuffer(tile1, tileinds)
imfilter!(r1, tileb1, A, k1, border, tileinds)
_imfilter_tiled_swap!(r1, out, kt, border, (tileb1, tile2))
end
out
end
# The first of the pair in `tmp` has the current data. We also make
# the second a plain array so there's no doubt about who's holding the
# proper indices.
function _imfilter_tiled_swap!(r, out, kernel::Tuple{Any,Any,Vararg{Any}}, border, tmp::Tuple{TileBuffer,Array})
tileb1, tile2 = tmp
k1, kt = kernel[1], tail(kernel)
parentinds = axes(tileb1)
tileinds = shrink(parentinds, k1)
tileb2 = TileBuffer(tile2, tileinds)
imfilter!(r, tileb2, tileb1, samedims(tileb2, k1), border, tileinds)
_imfilter_tiled_swap!(r, out, kt, border, (tileb2, parent(tileb1)))
end
# on the last call we write to `out` instead of one of the buffers
function _imfilter_tiled_swap!(r, out, kernel::Tuple{Any}, border, tmp::Tuple{TileBuffer,Array})
tileb1 = tmp[1]
k1 = kernel[1]
parentinds = axes(tileb1)
tileinds = shrink(parentinds, k1)
imfilter!(r, out, tileb1, samedims(out, k1), border, tileinds)
end
### FIR filtering
"""
imfilter!(::AbstractResource, imgfilt, img, kernel, NoPad(), [inds=axes(imgfilt)])
Filter an array `img` with kernel `kernel` by computing their
correlation, storing the result in `imgfilt`, defaulting to a finite-impulse
response (FIR) algorithm. Any necessary padding must have already been
supplied to `img`. If you want padding applied, instead call
imfilter!([r::AbstractResource,] imgfilt, img, kernel, border)
with a specific `border`, or use
imfilter!(imgfilt, img, kernel, [Algorithm.FIR()])
for default padding.
If `inds` is supplied, only the elements of `imgfilt` with indices in
the domain of `inds` will be calculated. This can be particularly
useful for "cascaded FIR filters" where you pad over a larger area and
then calculate the result over just the necessary/well-defined region
at each successive stage.
See also: [`imfilter`](@ref).
"""
function imfilter!(r::AbstractResource,
out::AbstractArray{S,N},
A::AbstractArray{T,N},
kern::NDimKernel{N},
border::NoPad,
inds::Indices{N}=axes(out)) where {S,T,N}
(isempty(A) || isempty(kern)) && return out
indso, indsA, indsk = axes(out), axes(A), axes(kern)
if iscopy(kern)
R = CartesianIndices(inds)
return copyto!(out, R, A, R)
end
for i = 1:N
# Check that inds is inbounds for out
indsi, indsoi, indsAi, indski = inds[i], indso[i], indsA[i], indsk[i]
if first(indsi) < first(indsoi) || last(indsi) > last(indsoi)
throw(DimensionMismatch("output indices $indso disagrees with requested indices $inds"))
end
# Check that input A is big enough not to throw a BoundsError
if first(indsAi) > first(indsi) + first(indski) ||
last(indsA[i]) < last(indsi) + last(indski)
throw(DimensionMismatch("requested indices $inds and kernel indices $indsk do not agree with indices of padded input, $indsA"))
end
end
_imfilter_inbounds!(r, out, A, kern, border, inds)
end
function _imfilter_inbounds!(r::AbstractResource, out, A::AbstractArray, kern::ReshapedIIR, border::NoPad, inds)
indspre, ind, indspost = iterdims(inds, kern)
_imfilter_dim!(r, out, A, kern.data, indspre, ind, indspost, border[])
end
function _imfilter_inbounds!(r::AbstractResource, out, A::AbstractArray, kern, border::NoPad, inds)
indsk = axes(kern)
R, Rk = CartesianIndices(inds), CartesianIndices(indsk)
if isempty(R) || isempty(Rk)
return out
end
p = accumfilter(A[first(R)+first(Rk)], first(kern))
z = zero(typeof(p+p))
__imfilter_inbounds!(r, out, A, kern, border, R, z)
end
function __imfilter_inbounds!(r, out, A, kern, border, R, z)
Rk = CartesianIndices(axes(kern))
for I in safetail(R), i in safehead(R)
tmp = z