This repository has been archived by the owner on Mar 3, 2021. It is now read-only.
forked from quatrope/astroalign
-
Notifications
You must be signed in to change notification settings - Fork 1
/
astroalign.py
482 lines (390 loc) · 18.1 KB
/
astroalign.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
# MIT License
# Copyright (c) 2016-2019 Martin Beroiz
# 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.
"""
ASTROALIGN is a simple package that will try to align two stellar astronomical
images, especially when there is no WCS information available.
It does so by finding similar 3-point asterisms (triangles) in both images and
deducing the affine transformation between them.
General registration routines try to match feature points, using corner
detection routines to make the point correspondence.
These generally fail for stellar astronomical images, since stars have very
little stable structure and so, in general, indistinguishable from each other.
Asterism matching is more robust, and closer to the human way of matching
stellar images.
Astroalign can match images of very different field of view, point-spread
functions, seeing and atmospheric conditions.
(c) Martin Beroiz
"""
__version__ = '2.0.1'
__all__ = [
"MAX_CONTROL_POINTS",
"MIN_MATCHES_FRACTION",
"MaxIterError",
"NUM_NEAREST_NEIGHBORS",
"PIXEL_TOL",
"apply_transform",
"estimate_transform",
"find_transform",
"matrix_transform",
"register"]
import numpy as _np
from skimage.transform import estimate_transform
from skimage.transform import matrix_transform # noqa
MAX_CONTROL_POINTS = 50
"""The maximum control points (stars) to use to build the invariants.
Default: 50"""
PIXEL_TOL = 2
"""The pixel distance tolerance to assume two invariant points are the same.
Default: 2"""
MIN_MATCHES_FRACTION = 0.8
"""The minimum fraction of triangle matches to accept a transformation.
If the minimum fraction yields more than 10 triangles, 10 is used instead.
Default: 0.8
"""
NUM_NEAREST_NEIGHBORS = 5
"""
The number of nearest neighbors of a given star (including itself) to construct
the triangle invariants.
Default: 5
"""
def _invariantfeatures(x1, x2, x3):
"Given 3 points x1, x2, x3, return the invariant features for the set."
sides = _np.sort([_np.linalg.norm(x1 - x2), _np.linalg.norm(x2 - x3),
_np.linalg.norm(x1 - x3)])
return [sides[2] / sides[1], sides[1] / sides[0]]
def _arrangetriplet(sources, vertex_indices):
"""Return vertex_indices ordered in an (a, b, c) form where:
a is the vertex defined by L1 & L2
b is the vertex defined by L2 & L3
c is the vertex defined by L3 & L1
and L1 < L2 < L3 are the sides of the triangle defined by vertex_indices."""
ind1, ind2, ind3 = vertex_indices
x1, x2, x3 = sources[vertex_indices]
side_ind = _np.array([(ind1, ind2), (ind2, ind3), (ind3, ind1)])
side_lengths = list(map(_np.linalg.norm, (x1 - x2, x2 - x3, x3 - x1)))
l1_ind, l2_ind, l3_ind = _np.argsort(side_lengths)
# the most common vertex in the list of vertices for two sides is the
# point at which they meet.
from collections import Counter
count = Counter(side_ind[[l1_ind, l2_ind]].flatten())
a = count.most_common(1)[0][0]
count = Counter(side_ind[[l2_ind, l3_ind]].flatten())
b = count.most_common(1)[0][0]
count = Counter(side_ind[[l3_ind, l1_ind]].flatten())
c = count.most_common(1)[0][0]
return _np.array([a, b, c])
def _generate_invariants(sources):
"""Return an array of (unique) invariants derived from the array `sources`.
Return an array of the indices of `sources` that correspond to each invariant,
arranged as described in _arrangetriplet.
"""
from scipy.spatial import KDTree
from itertools import combinations
from functools import partial
arrange = partial(_arrangetriplet, sources=sources)
inv = []
triang_vrtx = []
coordtree = KDTree(sources)
# The number of nearest neighbors to request (to work with few sources)
knn = min(len(sources), NUM_NEAREST_NEIGHBORS)
for asrc in sources:
__, indx = coordtree.query(asrc, knn)
# Generate all possible triangles with the 5 indx provided, and store
# them with the order (a, b, c) defined in _arrangetriplet
all_asterism_triang = [arrange(vertex_indices=list(cmb))
for cmb in combinations(indx, 3)]
triang_vrtx.extend(all_asterism_triang)
inv.extend([_invariantfeatures(*sources[triplet])
for triplet in all_asterism_triang])
# Remove here all possible duplicate triangles
uniq_ind = [pos for (pos, elem) in enumerate(inv)
if elem not in inv[pos + 1:]]
inv_uniq = _np.array(inv)[uniq_ind]
triang_vrtx_uniq = _np.array(triang_vrtx)[uniq_ind]
return inv_uniq, triang_vrtx_uniq
class _MatchTransform:
def __init__(self, source, target):
self.source = source
self.target = target
def fit(self, data):
"""
Return the best 2D similarity transform from the points given in data.
data: N sets of similar corresponding triangles.
3 indices for a triangle in ref
and the 3 indices for the corresponding triangle in target;
arranged in a (N, 3, 2) array.
"""
d1, d2, d3 = data.shape
s, d = data.reshape(d1 * d2, d3).T
approx_t = estimate_transform('similarity',
self.source[s], self.target[d])
return approx_t
def get_error(self, data, approx_t):
d1, d2, d3 = data.shape
s, d = data.reshape(d1 * d2, d3).T
resid = approx_t.residuals(self.source[s], self.target[d])\
.reshape(d1, d2)
error = resid.max(axis=1)
return error
def find_transform(source, target):
"""Estimate the transform between ``source`` and ``target``.
Return a SimilarityTransform object ``T`` that maps pixel x, y indices from
the source image s = (x, y) into the target (destination) image t = (x, y).
T contains parameters of the tranformation: ``T.rotation``,
``T.translation``, ``T.scale``, ``T.params``.
Args:
source (array-like): Either a numpy array of the source image to be
transformed or an interable of (x, y) coordinates of the target
control points.
target (array-like): Either a numpy array of the target (destination)
image or an interable of (x, y) coordinates of the target
control points.
Returns:
The transformation object and a tuple of corresponding star positions
in source and target.::
T, (source_pos_array, target_pos_array)
Raises:
TypeError: If input type of ``source`` or ``target`` is not supported.
Exception: If it cannot find more than 3 stars on any input.
"""
from scipy.spatial import KDTree
try:
if len(source[0]) == 2:
# Assume it's a list of (x, y) pairs
source_controlp = _np.array(source)[:MAX_CONTROL_POINTS]
else:
# Assume it's a 2D image
source_controlp = _find_sources(source)[:MAX_CONTROL_POINTS]
except:
raise TypeError('Input type for source not supported.')
try:
if len(target[0]) == 2:
# Assume it's a list of (x, y) pairs
target_controlp = _np.array(target)[:MAX_CONTROL_POINTS]
else:
# Assume it's a 2D image
target_controlp = _find_sources(target)[:MAX_CONTROL_POINTS]
except:
raise TypeError('Input type for target not supported.')
# Check for low number of reference points
if len(source_controlp) < 3:
raise Exception("Reference stars in source image are less than the "
"minimum value (3).")
if len(target_controlp) < 3:
raise Exception("Reference stars in target image are less than the "
"minimum value (3).")
source_invariants, source_asterisms = _generate_invariants(source_controlp)
source_invariant_tree = KDTree(source_invariants)
target_invariants, target_asterisms = _generate_invariants(target_controlp)
target_invariant_tree = KDTree(target_invariants)
# r = 0.1 is the maximum search distance, 0.1 is an empirical value that
# returns about the same number of matches than inputs
# matches_list is a list of lists such that for each element
# source_invariant_tree.data[i], matches_list[i] is a list of the indices
# of its neighbors in target_invariant_tree.data
matches_list = \
source_invariant_tree.query_ball_tree(target_invariant_tree, r=0.1)
# matches unravels the previous list of matches into pairs of source and
# target control point matches.
# matches is a (N, 3, 2) array. N sets of similar corresponding triangles.
# 3 indices for a triangle in ref
# and the 3 indices for the corresponding triangle in target;
matches = []
# t1 is an asterism in source, t2 in target
for t1, t2_list in zip(source_asterisms, matches_list):
for t2 in target_asterisms[t2_list]:
matches.append(list(zip(t1, t2)))
matches = _np.array(matches)
inv_model = _MatchTransform(source_controlp, target_controlp)
n_invariants = len(matches)
max_iter = n_invariants
# Set the minimum matches to be between 1 and 10 asterisms
min_matches = max(1, min(10, int(n_invariants * MIN_MATCHES_FRACTION)))
if (len(source_controlp) == 3 or len(target_controlp) == 3)\
and len(matches) == 1:
best_t = inv_model.fit(matches)
inlier_ind = _np.arange(len(matches)) # All of the indices
else:
best_t, inlier_ind = _ransac(
matches, inv_model, 1, max_iter, PIXEL_TOL, min_matches
)
triangle_inliers = matches[inlier_ind]
d1, d2, d3 = triangle_inliers.shape
inl_arr = triangle_inliers.reshape(d1 * d2, d3)
inl_unique = set(tuple(pair) for pair in inl_arr)
inl_arr_unique = _np.array(list(list(apair) for apair in inl_unique))
s, d = inl_arr_unique.T
return best_t, (source_controlp[s], target_controlp[d])
def apply_transform(transform, source, target,
fill_value=None, propagate_mask=False):
"""Applies the transformation ``transform`` to ``source``.
The output image will have the same shape as ``target``.
Args:
transform: A scikit-image ``SimilarityTransform`` object.
source (numpy array): A 2D numpy array of the source image to be
transformed.
target (numpy array): A 2D numpy array of the target image. Only used
to set the output image shape.
fill_value (float): A value to fill in the areas of aligned_image
where footprint == True.
propagate_mask (bool): Wether to propagate the mask in source.mask
onto footprint.
Return:
A tuple (aligned_image, footprint).
aligned_image is a numpy 2D array of the transformed source
footprint is a mask 2D array with True on the regions
with no pixel information.
"""
from skimage.transform import warp
if hasattr(source, 'data') and isinstance(source.data, _np.ndarray):
source_data = source.data
else:
source_data = source
if hasattr(target, 'data') and isinstance(target.data, _np.ndarray):
target_data = target.data
else:
target_data = target
aligned_image = warp(source_data, inverse_map=transform.inverse,
output_shape=target_data.shape, order=3, mode='constant',
cval=_np.median(source_data), clip=False,
preserve_range=True)
footprint = warp(_np.zeros(source_data.shape, dtype='float32'),
inverse_map=transform.inverse,
output_shape=target_data.shape,
cval=1.0)
footprint = footprint > 0.4
if hasattr(source, 'mask') and propagate_mask:
source_mask = _np.array(source.mask)
if source_mask.shape == source_data.shape:
source_mask_rot = warp(source_mask.astype('float32'),
inverse_map=transform.inverse,
output_shape=target_data.shape,
cval=1.0)
source_mask_rot = source_mask_rot > 0.4
footprint = footprint | source_mask_rot
if fill_value is not None:
aligned_image[footprint] = fill_value
return aligned_image, footprint
def register(source, target, fill_value=None, propagate_mask=False):
"""Transform ``source`` to coincide pixel to pixel with ``target``.
Args:
source (numpy array): A 2D numpy array of the source image to be
transformed.
target (numpy array): A 2D numpy array of the target image. Only used
to set the output image shape.
fill_value (float): A value to fill in the areas of aligned_image
where footprint == True.
propagate_mask (bool): Wether to propagate the mask in source.mask
onto footprint.
Return:
A tuple (aligned_image, footprint).
aligned_image is a numpy 2D array of the transformed source
footprint is a mask 2D array with True on the regions
with no pixel information.
"""
t, __ = find_transform(source=source, target=target)
aligned_image, footprint = apply_transform(t, source, target,
fill_value, propagate_mask,
)
return aligned_image, footprint
def _find_sources(img):
"Return sources (x, y) sorted by brightness."
import sep
if isinstance(img, _np.ma.MaskedArray):
image = img.filled(fill_value=_np.median(img)).astype('float32')
else:
image = img.astype('float32')
bkg = sep.Background(image)
thresh = 3. * bkg.globalrms
sources = sep.extract(image - bkg.back(), thresh)
sources.sort(order='flux')
return _np.array([[asrc['x'], asrc['y']] for asrc in sources[::-1]])
# Copyright (c) 2004-2007, Andrew D. Straw. All rights reserved.
# Redistribution and use in source and binary forms, with or without
# modification, are permitted provided that the following conditions are
# met:
# * Redistributions of source code must retain the above copyright
# notice, this list of conditions and the following disclaimer.
# * Redistributions in binary form must reproduce the above
# copyright notice, this list of conditions and the following
# disclaimer in the documentation and/or other materials provided
# with the distribution.
# * Neither the name of the Andrew D. Straw nor the names of its
# contributors may be used to endorse or promote products derived
# from this software without specific prior written permission.
# THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
# "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
# LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
# a PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
# OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
# SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
# LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
# DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
# THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
# (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
# OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
#
#
# Modified by Martin Beroiz
class MaxIterError(Exception):
pass
def _ransac(data, model, min_data_points, max_iter, thresh, min_matches):
"""fit model parameters to data using the RANSAC algorithm
This implementation written from pseudocode found at
http://en.wikipedia.org/w/index.php?title=RANSAC&oldid=116358182
Given:
data: a set of data points
model: a model that can be fitted to data points
min_data_points: the minimum number of data values required to fit the
model
max_iter: the maximum number of iterations allowed in the algorithm
thresh: a threshold value to determine when a data point fits a model
min_matches: the min number of matches required to assert that a model
fits well to data
Return:
bestfit: model parameters which best fit the data (or nil if no good model
is found)
"""
iterations = 0
bestfit = None
best_inlier_idxs = None
n_data = data.shape[0]
n = min_data_points
all_idxs = _np.arange(n_data)
while iterations < max_iter:
# Partition indices into two random subsets
_np.random.shuffle(all_idxs)
maybe_idxs, test_idxs = all_idxs[:n], all_idxs[n:]
maybeinliers = data[maybe_idxs, :]
test_points = data[test_idxs, :]
maybemodel = model.fit(maybeinliers)
test_err = model.get_error(test_points, maybemodel)
# select indices of rows with accepted points
also_idxs = test_idxs[test_err < thresh]
alsoinliers = data[also_idxs, :]
if len(alsoinliers) >= min_matches:
betterdata = _np.concatenate((maybeinliers, alsoinliers))
bestfit = model.fit(betterdata)
best_inlier_idxs = _np.concatenate((maybe_idxs, also_idxs))
break
iterations += 1
if bestfit is None:
raise MaxIterError(
"Max iterations exceeded while trying to find "
"acceptable transformation.")
return bestfit, best_inlier_idxs