-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathsolver_utils.py
286 lines (253 loc) · 8.25 KB
/
solver_utils.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
import numpy as np
from modopt.opt.algorithms import POGM, ForwardBackward, Condat
from modopt.opt.linear import Identity
from mri.operators.gradient.gradient import GradAnalysis, GradSynthesis
OPTIMIZERS = {
"pogm": "synthesis",
"fista": "analysis",
"condat-vu": "analysis",
None: None,
}
def get_grad_op(fourier_op, grad_formulation, linear_op=None, verbose=False, **kwargs):
"""Create gradient operator specific to the problem."""
if grad_formulation == "analysis":
return GradAnalysis(fourier_op=fourier_op, verbose=verbose, **kwargs)
if grad_formulation == "synthesis":
return GradSynthesis(
linear_op=linear_op,
fourier_op=fourier_op,
verbose=verbose,
**kwargs,
)
def initialize_opt(
opt_name,
grad_op,
linear_op,
prox_op,
x_init=None,
synthesis_init=False,
opt_kwargs=None,
metric_kwargs=None,
):
"""
Initialize an Optimizer with the suitable parameters.
Parameters:
----------
grad_op: OperatorBase
Gradient Operator for the data consistency
x_init: ndarray, default None
Initial value for the reconstruction. If None use a zero Array.
synthesis_init: bool, default False
Is the initial_value in the image space of the space_linear operator ?
opt_kwargs: dict, default None
Extra kwargs for the initialisation of Optimizer
metric_kwargs: dict, default None
Extra kwargs for the metric api of ModOpt
Returns:
-------
An Optimizer Instance to perform the reconstruction with.
See Also:
--------
Modopt.opt.algorithms
"""
if x_init is None:
x_init = np.squeeze(
np.zeros(
(grad_op.fourier_op.n_coils, *grad_op.fourier_op.shape),
dtype="complex64",
)
)
if not synthesis_init and hasattr(grad_op, "linear_op"):
alpha_init = grad_op.linear_op.op(x_init)
elif synthesis_init and not hasattr(grad_op, "linear_op"):
x_init = linear_op.adj_op(x_init)
elif not synthesis_init and hasattr(grad_op, "linear_op"):
alpha_init = x_init
opt_kwargs = opt_kwargs or dict()
metric_kwargs = metric_kwargs or dict()
beta = grad_op.inv_spec_rad
if opt_name == "pogm":
opt = POGM(
u=alpha_init,
x=alpha_init,
y=alpha_init,
z=alpha_init,
grad=grad_op,
prox=prox_op,
linear=linear_op,
beta_param=beta,
sigma_bar=opt_kwargs.pop("sigma_bar", 0.96),
auto_iterate=opt_kwargs.pop("auto_iterate", False),
**opt_kwargs,
**metric_kwargs,
)
elif opt_name == "fista":
opt = ForwardBackward(
x=x_init,
grad=grad_op,
prox=prox_op,
linear=linear_op,
beta_param=beta,
lambda_param=opt_kwargs.pop("lambda_param", 1.0),
auto_iterate=opt_kwargs.pop("auto_iterate", False),
**opt_kwargs,
**metric_kwargs,
)
elif opt_name == "condat-vu":
y_init = linear_op.op(x_init)
opt = Condat(
x=x_init,
y=y_init,
grad=grad_op,
prox=Identity(),
prox_dual=prox_op,
linear=linear_op,
**opt_kwargs,
**metric_kwargs,
)
else:
raise ValueError(f"Optimizer {opt_name} not implemented")
return opt
from modopt.opt.linear import LinearParent
import pywt
from joblib import Parallel, delayed, cpu_count
import numpy as np
class WaveletTransform(LinearParent):
"""
2D and 3D wavelet transform class.
This is a light wrapper around PyWavelet, with multicoil support.
Parameters
----------
wavelet_name: str
the wavelet name to be used during the decomposition.
shape: tuple[int,...]
Shape of the input data. The shape should be a tuple of length 2 or 3.
It should not contains coils or batch dimension.
nb_scales: int, default 4
the number of scales in the decomposition.
n_coils: int, default 1
the number of coils for multichannel reconstruction
n_jobs: int, default 1
the number of cores to use for multichannel.
backend: str, default "threading"
the backend to use for parallel multichannel linear operation.
verbose: int, default 0
the verbosity level.
Attributes
----------
nb_scale: int
number of scale decomposed in wavelet space.
n_jobs: int
number of jobs for parallel computation
n_coils: int
number of coils use f
backend: str
Backend use for parallel computation
verbose: int
Verbosity level
"""
def __init__(
self,
wavelet_name,
shape,
level=4,
n_coils=1,
n_jobs=1,
decimated=True,
backend="threading",
mode="symmetric",
):
if wavelet_name not in pywt.wavelist(kind="all"):
raise ValueError(
"Invalid wavelet name. Check ``pywt.waveletlist(kind='all')``"
)
self.wavelet = wavelet_name
if isinstance(shape, int):
shape = (shape,)
self.shape = shape
self.n_jobs = n_jobs
self.mode = mode
self.level = level
if not decimated:
raise NotImplementedError(
"Undecimated Wavelet Transform is not implemented yet."
)
ca, *cds = pywt.wavedecn_shapes(
self.shape, wavelet=self.wavelet, mode=self.mode, level=self.level
)
self.coeffs_shape = [ca] + [s for cd in cds for s in cd.values()]
if len(shape) > 1:
self.dwt = pywt.wavedecn
self.idwt = pywt.waverecn
self._pywt_fun = "wavedecn"
else:
self.dwt = pywt.wavedec
self.idwt = pywt.waverec
self._pywt_fun = "wavedec"
self.n_coils = n_coils
if self.n_coils == 1 and self.n_jobs != 1:
print("Making n_jobs = 1 for WaveletN as n_coils = 1")
self.n_jobs = 1
self.backend = backend
n_proc = self.n_jobs
if n_proc < 0:
n_proc = cpu_count() + self.n_jobs + 1
def op(self, data):
"""Define the wavelet operator.
This method returns the input data convolved with the wavelet filter.
Parameters
----------
data: ndarray or Image
input 2D data array.
Returns
-------
coeffs: ndarray
the wavelet coefficients.
"""
if self.n_coils > 1:
coeffs, self.coeffs_slices, self.raw_coeffs_shape = zip(
*Parallel(
n_jobs=self.n_jobs, backend=self.backend, verbose=self.verbose
)(delayed(self._op)(data[i]) for i in np.arange(self.n_coils))
)
coeffs = np.asarray(coeffs)
else:
coeffs, self.coeffs_slices, self.raw_coeffs_shape = self._op(data)
return coeffs
def _op(self, data):
"""Single coil wavelet transform."""
return pywt.ravel_coeffs(
self.dwt(data, mode=self.mode, level=self.level, wavelet=self.wavelet)
)
def adj_op(self, coeffs):
"""Define the wavelet adjoint operator.
This method returns the reconstructed image.
Parameters
----------
coeffs: ndarray
the wavelet coefficients.
Returns
-------
data: ndarray
the reconstructed data.
"""
if self.n_coils > 1:
images = Parallel(
n_jobs=self.n_jobs, backend=self.backend, verbose=self.verbose
)(
delayed(self._adj_op)(coeffs[i], self.coeffs_shape[i])
for i in np.arange(self.n_coils)
)
images = np.asarray(images)
else:
images = self._adj_op(coeffs)
return images
def _adj_op(self, coeffs):
"""Single coil inverse wavelet transform."""
return self.idwt(
pywt.unravel_coeffs(
coeffs, self.coeffs_slices, self.raw_coeffs_shape, self._pywt_fun
),
wavelet=self.wavelet,
mode=self.mode,
)