/
peakfinders2D.py
320 lines (273 loc) · 10.9 KB
/
peakfinders2D.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
# -*- coding: utf-8 -*-
# Copyright 2017-2018 The pyXem developers
#
# This file is part of pyXem.
#
# pyXem is free software: you can redistribute it and/or modify
# it under the terms of the GNU General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# pyXem is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU General Public License for more details.
#
# You should have received a copy of the GNU General Public License
# along with pyXem. If not, see <http://www.gnu.org/licenses/>.
import numpy as np
import scipy.ndimage as ndi
NO_PEAKS = np.array([[[np.nan, np.nan]]])
def clean_peaks(peaks):
if len(peaks) == 0:
return NO_PEAKS
else:
return peaks
def find_peaks_zaefferer(z, grad_threshold=0.1, window_size=40,
distance_cutoff=50.):
"""Method to locate positive peaks in an image based on gradient
thresholding and subsequent refinement within masked regions.
Parameters
----------
z : numpy.ndarray
Matrix of image intensities.
grad_threshold : float
The minimum gradient required to begin a peak search.
window_size : int
The size of the square window within which a peak search is
conducted. If odd, will round down to even.
distance_cutoff : float
The maximum distance a peak may be from the initial
high-gradient point.
Returns
-------
peaks : numpy.ndarray
(n_peaks, 2)
Peak pixel coordinates.
Notes
-----
Implemented as described in Zaefferer "New developments of computer-aided
crystallographic analysis in transmission electron microscopy" J. Ap. Cryst.
This version by Ben Martineau (2016)
"""
def box(x, y, window_size, x_max, y_max):
"""Produces a list of coordinates in the box about (x, y)."""
a = int(window_size / 2)
x_min = max(0, x - a)
x_max = min(x_max, x + a)
y_min = max(0, y - a)
y_max = min(y_max, y + a)
return np.array(
np.meshgrid(range(x_min, x_max), range(y_min, y_max))).reshape(
2, -1).T
def get_max(image, box):
"""Finds the coordinates of the maximum of 'image' in 'box'."""
vals = image[box[:, 0], box[:, 1]]
max_position = box[np.argmax(vals)]
return max_position
def distance(x, y):
"""Calculates the distance between two points."""
v = x - y
return np.sqrt(np.sum(np.square(v)))
def gradient(image):
"""Calculates the square of the 2-d partial gradient.
Parameters
----------
image : numpy.ndarray
Returns
-------
numpy.ndarray
"""
gradient_of_image = np.gradient(image)
gradient_of_image = gradient_of_image[0] ** 2 + gradient_of_image[
1] ** 2
return gradient_of_image
# Generate an ordered list of matrix coordinates.
z = z / np.max(z)
coordinates = np.indices(z.data.shape).reshape(2, -1).T
# Calculate the gradient at every point.
image_gradient = gradient(z)
# Boolean matrix of high-gradient points.
gradient_is_above_threshold = image_gradient >= grad_threshold
peaks = []
for coordinate in coordinates[gradient_is_above_threshold.flatten()]:
# Iterate over coordinates where the gradient is high enough.
b = box(coordinate[0], coordinate[1], window_size, z.shape[0],
z.shape[1])
p_old = np.array([0, 0])
p_new = get_max(z, b)
while np.all(p_old != p_new):
p_old = p_new
b = box(p_old[0], p_old[1], window_size, z.shape[0], z.shape[1])
p_new = get_max(z, b)
if distance(coordinate, p_new) > distance_cutoff:
break
peaks.append(tuple(p_new))
peaks = np.array([np.array(p) for p in set(peaks)])
return clean_peaks(peaks)
def find_peaks_stat(z, alpha=1., window_radius=10, convergence_ratio=0.05):
"""Locate positive peaks in an image based on statistical refinement and
difference with respect to mean intensity.
Parameters
----------
z : numpy.ndarray
Array of image intensities.
alpha : float
Only maxima above `alpha * sigma` are found, where `sigma` is the
local, rolling standard deviation of the image.
window_radius : int
The pixel radius of the circular window for the calculation of the
rolling mean and standard deviation.
convergence_ratio : float
The algorithm will stop finding peaks when the proportion of new peaks
being found is less than `convergence_ratio`.
Returns
-------
numpy.ndarray
(n_peaks, 2)
Array of peak coordinates.
Notes
-----
Implemented as described in the PhD thesis of Thomas White (2009) the
algorithm was developed by Gordon Ball during a summer project in
Cambridge.
This version by Ben Martineau (2016), with minor modifications to the
original where methods were ambiguous or unclear.
"""
from scipy.ndimage.filters import generic_filter
from scipy.ndimage.filters import uniform_filter
from sklearn.cluster import DBSCAN
def normalize(image):
"""Scales the image to intensities between 0 and 1."""
return image / np.max(image)
def _local_stat(image, radius, func):
"""Calculates rolling method 'func' over a circular kernel."""
x, y = np.ogrid[-radius:radius + 1, -radius:radius + 1]
kernel = np.hypot(x, y) < radius
stat = generic_filter(image, func, footprint=kernel)
return stat
def local_mean(image, radius):
"""Calculates rolling mean over a circular kernel."""
return _local_stat(image, radius, np.mean)
def local_std(image, radius):
"""Calculates rolling standard deviation over a circular kernel."""
return _local_stat(image, radius, np.std)
def single_pixel_desensitize(image):
"""Reduces single-pixel anomalies by nearest-neighbor smoothing."""
kernel = np.array([[0.5, 1, 0.5], [1, 1, 1], [0.5, 1, 0.5]])
smoothed_image = generic_filter(image, np.mean, footprint=kernel)
return smoothed_image
def stat_binarise(image):
"""Peaks more than one standard deviation from the mean set to one."""
image_rolling_mean = local_mean(image, window_radius)
image_rolling_std = local_std(image, window_radius)
image = single_pixel_desensitize(image)
binarised_image = np.zeros(image.shape)
stat_mask = image > (image_rolling_mean + alpha * image_rolling_std)
binarised_image[stat_mask] = 1
return binarised_image
def smooth(image):
"""Image convolved twice using a uniform 3x3 kernel."""
image = uniform_filter(image, size=3)
image = uniform_filter(image, size=3)
return image
def half_binarise(image):
"""Image binarised about values of one-half intensity."""
binarised_image = np.where(image > 0.5, 1, 0)
return binarised_image
def separate_peaks(binarised_image):
"""Identify adjacent 'on' coordinates via DBSCAN."""
bi = binarised_image.astype('bool')
coordinates = np.indices(bi.shape).reshape(2, -1).T[
bi.flatten()]
db = DBSCAN(2, 3)
peaks = []
labeled_points = db.fit_predict(coordinates)
for peak_label in list(set(labeled_points)):
peaks.append(coordinates[labeled_points == peak_label])
return peaks
def _peak_find_once(image):
"""Smooth, binarise, and find peaks according to main algorithm."""
image = smooth(image)
image = half_binarise(image)
peaks = separate_peaks(image)
return image, peaks
def stat_peak_finder(image):
"""Find peaks in image. Algorithm stages in comments."""
image = normalize(image) # 1
image = stat_binarise(image) # 2, 3
n_peaks = np.infty # Initial number of peaks
image, peaks = _peak_find_once(image) # 4-6
m_peaks = len(peaks) # Actual number of peaks
"""
#XXX
while (n_peaks - m_peaks) / n_peaks > convergence_ratio: # 8
n_peaks = m_peaks
image, peaks = _peak_find_once(image)
m_peaks = len(peaks)
"""
peak_centers = np.array(
[np.mean(peak, axis=0) for peak in peaks]) # 7
return peak_centers
return clean_peaks(stat_peak_finder(z))
def find_peaks_dog(z, min_sigma=1., max_sigma=50., sigma_ratio=1.6,
threshold=0.2, overlap=0.5):
"""
Finds peaks via the difference of Gaussian Matrices method from
`scikit-image`.
Parameters
----------
z : numpy.ndarray
2-d array of intensities
float min_sigma, max_sigma, sigma_ratio, threshold, overlap
Additional parameters to be passed to the algorithm. See `blob_dog`
documentation for details:
http://scikit-image.org/docs/dev/api/skimage.feature.html#blob-dog
Returns
-------
numpy.ndarray
Array of peak coordinates of shape `(n_peaks, 2)`
Notes
-----
While highly effective at finding even very faint peaks, this method is
sensitive to fluctuations in intensity near the edges of the image.
"""
from skimage.feature import blob_dog
z = z / np.max(z)
blobs = blob_dog(z, min_sigma=min_sigma, max_sigma=max_sigma,
sigma_ratio=sigma_ratio, threshold=threshold,
overlap=overlap)
centers = blobs[:, :2]
clean_centers = []
for center in centers:
if len(np.intersect1d(center, (0, 1) + z.shape + tuple(
c - 1 for c in z.shape))) > 0:
continue
clean_centers.append(center)
return np.array(clean_centers)
def find_peaks_log(z, min_sigma=1., max_sigma=50., num_sigma=10.,
threshold=0.2, overlap=0.5, log_scale=False):
"""
Finds peaks via the Laplacian of Gaussian Matrices method from
`scikit-image`.
Parameters
----------
z : numpy.ndarray
Array of image intensities.
float min_sigma, max_sigma, num_sigma, threshold, overlap, log_scale
Additional parameters to be passed to the algorithm. See
`blob_log` documentation for details:
http://scikit-image.org/docs/dev/api/skimage.feature.html#blob-log
Returns
-------
numpy.ndarray
(n_peaks, 2)
Array of peak coordinates.
"""
from skimage.feature import blob_log
z = z / np.max(z)
blobs = blob_log(z, min_sigma=min_sigma, max_sigma=max_sigma,
num_sigma=num_sigma, threshold=threshold, overlap=overlap,
log_scale=log_scale)
centers = blobs[:, :2]
return centers