-
Notifications
You must be signed in to change notification settings - Fork 21
/
mean.py
323 lines (258 loc) · 11.3 KB
/
mean.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
import logging
from functools import partial
import numpy as np
from scipy.linalg import norm
from scipy.sparse.linalg import LinearOperator, cg
from aspire import config
from aspire.basis import Coef
from aspire.nufft import anufft
from aspire.numeric import fft
from aspire.operators import evaluate_src_filters_on_grid
from aspire.reconstruction import Estimator, FourierKernel, FourierKernelMatrix
from aspire.volume import Volume, rotated_grids
logger = logging.getLogger(__name__)
class WeightedVolumesEstimator(Estimator):
def __init__(self, weights, *args, **kwargs):
"""
Weighted mean volume estimation.
This class holds the `FourierKernelMatrix`, stored as a r x r
matrix of volumes. The problem being solved here is the
minimization given by eq. (14) in the paper, rewritten as the
normal equations in eq. (20) and more compactly in eq. (23).
Convolution with each of these kernels is equivalent to
performing a projection/backprojection on a volume, with the
appropriate amplitude modifiers and CTF, and also a weighting
term; the r^2 volumes are each of pairwise products between
the weighting vectors given by the columns of wts.
Note that this is a non-centered Fourier transform, so the
zero frequency is found at index 0.
Formulas and conventions used for Volume estimation are
described in:
“Structural Variability from Noisy Tomographic Projections.”
Andén, Joakim, and Amit Singer.
SIAM journal on imaging sciences vol. 11,2 (2018): 1441-1492.
doi:10.1137/17M1153509
"Cryo-EM reconstruction of continuous heterogeneity by Laplacian spectral volumes"
Amit Moscovich, Amit Halevi, Joakim Andén and Amit Singer
Inverse Problems, Volume 36, Number 2, 2020 IOP Publishing Ltd
Special Issue on Cryo-Electron Microscopy and Inverse Problems
https://doi.org/10.1088/1361-6420/ab4f55
:param weights: Matrix of weights, n x r.
"""
self.weights = weights
self.r = self.weights.shape[1]
super().__init__(*args, **kwargs)
assert self.src.n == self.weights.shape[0]
def __getattr__(self, name):
if name == "precond_kernel":
if self.preconditioner == "circulant":
self.precond_kernel = FourierKernelMatrix(
1.0 / self.kernel.circularize()
)
else:
if self.preconditioner.lower() not in (None, "none"):
logger.warning(
f"Preconditioner {self.preconditioner} is not implemented, resetting to default of None."
)
self.precond_kernel = None
return self.precond_kernel
else:
# raise AttributeError(name)
return super().__getattr__(name)
def compute_kernel(self):
"""
Compute and return FourierKernelMatrix instance.
"""
return FourierKernelMatrix(self._compute_kernel())
def _compute_kernel(self):
"""
:return: r x r matrix of volumes, shaped (r, r, 2L, 2L, 2L).
"""
_2L = 2 * self.src.L
# Note, because we're iteratively summing it is critical we zero this array.
kernel = np.zeros((self.r, self.r, _2L, _2L, _2L), dtype=self.dtype)
# Handle symmetry boosting.
sym_rots = np.eye(3, dtype=self.dtype)[None]
if self.boost:
sym_rots = self.src.symmetry_group.matrices
for i in range(0, self.src.n, self.batch_size):
_range = np.arange(i, min(self.src.n, i + self.batch_size), dtype=int)
sq_filters_f = evaluate_src_filters_on_grid(self.src, _range) ** 2
amplitudes_sq = (self.src.amplitudes[_range] ** 2).astype(
self.dtype, copy=False
)
for k in range(self.r):
for j in range(k + 1):
weights = sq_filters_f * amplitudes_sq
if self.src.L % 2 == 0:
weights[0, :, :] = 0
weights[:, 0, :] = 0
weights *= (
self.weights[_range, j] * self.weights[_range, k]
).reshape(1, 1, len(_range))
weights = np.transpose(weights, (2, 0, 1)).flatten()
# Apply boosting.
batch_kernel = np.zeros((_2L, _2L, _2L), dtype=self.dtype)
for sym_rot in sym_rots:
rotations = sym_rot @ self.src.rotations[_range]
pts_rot = rotated_grids(self.src.L, rotations)
pts_rot = pts_rot.reshape((3, -1))
batch_kernel += (
1
/ (self.r * self.src.L**4)
* anufft(weights, pts_rot[::-1], (_2L, _2L, _2L), real=True)
)
kernel[k, j] += batch_kernel
# r x r symmetric
# accumulate batch entries of kernel[k,j] to kernel[j,k]
if j != k:
kernel[j, k] += batch_kernel
kermat_f = np.zeros((self.r, self.r, _2L, _2L, _2L))
logger.info("Computing non-centered Fourier Transform Kernel Mat")
for k in range(self.r):
for j in range(self.r):
# Ensure symmetric kernel
kernel[k, j, 0, :, :] = 0
kernel[k, j, :, 0, :] = 0
kernel[k, j, :, :, 0] = 0
kernel[k, j] = fft.mdim_ifftshift(kernel[k, j], range(0, 3))
kernel_f = fft.fftn(kernel[k, j], axes=(0, 1, 2))
kernel_f = np.real(kernel_f)
kermat_f[k, j] = kernel_f
return kermat_f
def src_backward(self):
"""
Apply adjoint mapping to source
:return: The adjoint mapping applied to the images, averaged over the whole dataset and expressed
as coefficients of `basis`.
"""
# Handle symmetry boosting.
symmetry_group = None
sym_order = 1
if self.boost:
symmetry_group = self.src.symmetry_group
sym_order = len(symmetry_group.matrices)
# src_vols_wt_backward
vol_rhs = Volume(
np.zeros((self.r, self.src.L, self.src.L, self.src.L), dtype=self.dtype)
)
for i in range(0, self.src.n, self.batch_size):
for k in range(self.r):
im = self.src.images[i : i + self.batch_size]
batch_vol_rhs = self.src.im_backward(
im,
i,
self.weights[:, k],
symmetry_group=symmetry_group,
) / (self.src.n * sym_order)
vol_rhs[k] += batch_vol_rhs.astype(self.dtype)
res = np.sqrt(self.src.n * sym_order) * self.basis.evaluate_t(vol_rhs)
logger.info(f"Determined weighted adjoint mappings. Shape = {res.shape}")
return res
def conj_grad(self, b_coef, tol=1e-5, regularizer=0):
count = b_coef.shape[-1] # b_coef should be (r, basis.count)
kernel = self.kernel
if regularizer > 0:
kernel += regularizer
operator = LinearOperator(
(self.r * count, self.r * count),
matvec=partial(self.apply_kernel, kernel=kernel),
dtype=self.dtype,
)
if self.precond_kernel is None:
M = None
else:
precond_kernel = self.precond_kernel
if regularizer > 0:
precond_kernel += regularizer
M = LinearOperator(
(self.r * count, self.r * count),
matvec=partial(self.apply_kernel, kernel=precond_kernel),
dtype=self.dtype,
)
tol = tol or config.mean.cg_tol
target_residual = tol * norm(b_coef)
# callback setup
self.i = 0 # iteration counter
def cb(xk):
self.i += 1 # increment iteration count
logger.info(
f"[Iter {self.i}]: Delta {norm(b_coef - self.apply_kernel(xk))} (target {target_residual})"
)
# Do checkpoint at `checkpoint_iterations`,
_do_checkpoint = (
self.checkpoint_iterations
and (self.i % self.checkpoint_iterations) == 0
)
# or the last iteration when `maxiter` provided.
if self.maxiter:
_do_checkpoint |= self.i == (self.maxiter - 1)
# Optional checkpoint
if _do_checkpoint:
# Construct checkpoint path
path = f"{self.checkpoint_prefix}_iter{self.i:04d}.npy"
# Write out the current solution
np.save(path, xk.reshape(self.r, self.basis.count))
logger.info(f"Checkpoint saved to `{path}`")
x, info = cg(
operator,
b_coef.flatten(),
M=M,
callback=cb,
tol=tol,
atol=0,
maxiter=self.maxiter,
)
if info != 0:
raise RuntimeError("Unable to converge!")
return x.reshape(self.r, self.basis.count)
def apply_kernel(self, vol_coef, kernel=None):
"""
Applies the kernel represented by convolution
:param vol_coef: The volume to be convolved, stored in the basis coefficients.
:param kernel: a Kernel object. If None, the kernel for this Estimator is used.
:return: The result of evaluating `vol_coef` in the given basis, convolving with the kernel given by
kernel, and backprojecting into the basis.
"""
if kernel is None:
kernel = self.kernel
assert np.size(vol_coef) == self.r * self.basis.count
if vol_coef.ndim == 1:
vol_coef = vol_coef.reshape(self.r, self.basis.count)
vols_out = Volume(
np.zeros((self.r, self.src.L, self.src.L, self.src.L), dtype=self.dtype)
)
vol = Coef(self.basis, vol_coef).evaluate()
for k in range(self.r):
for j in range(self.r):
vols_out[k] = vols_out[k] + kernel.convolve_volume(vol[j], j, k)
# Note this is where we would add mask_gamma
vol_coef = self.basis.evaluate_t(vols_out)
return vol_coef
class MeanEstimator(WeightedVolumesEstimator):
"""
Special case of weighted mean volume estimate,
for a single volume.
"""
def __init__(self, src, **kwargs):
# Note, Handle boosting by adjusting weights based on symmetric order.
weights = np.ones((src.n, 1)) / np.sqrt(
src.n * len(src.symmetry_group.matrices)
)
super().__init__(weights, src, **kwargs)
def __getattr__(self, name):
"""
See `Estimator.__getattr__`.
"""
return super(WeightedVolumesEstimator, self).__getattr__(name)
def apply_kernel(self, *args, **kwargs):
"""
See `Estimator.apply_kernel`.
"""
return super(WeightedVolumesEstimator, self).apply_kernel(*args, **kwargs)
def compute_kernel(self):
"""
Compute and return `FourierKernel` instance.
"""
# Note for the r=1 we select and return a single kernel.
return FourierKernel(self._compute_kernel()[0][0])