-
Notifications
You must be signed in to change notification settings - Fork 0
/
kspace_det_edg.py
233 lines (185 loc) · 8.58 KB
/
kspace_det_edg.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
# -*- coding: utf-8 -*-
"""
k-space based details/edges detection in MRI images with optional k-space based
denoising and detail control (for Agilent FID data).
Created on Tue Nov 15 2022
Last modified on Thu Dec 1 2022
@author: Beata Wereszczyńska
"""
import nmrglue as ng
import numpy as np
import cv2
from skimage import morphology
import matplotlib.pyplot as plt
def import_kspace(path, number_of_slices, picked_slice):
"""
Imports k-space corresponding to the picked slice from Agilent FID data.
Input:
.fid folder location: path [str],
total number of slices in the MRI experiment: number_of_slices [int],
selected slice number: picked_slice [int],
"""
# import k-space data
echoes = ng.agilent.read(dir=path)[1]
kspace = echoes[picked_slice - 1 : echoes.shape[0] : number_of_slices, :] # downsampling to one slice
# return data
return kspace
def wght_msk_kspace(kspace, weight_power, contrast):
"""
k-space weighting and masking for denoising of MRI image without blurring or losing contrast,
as well as for brightening of the objects in the image with simultaneous noise reduction.
Input:
k-space data: kspace [array]
exponent in the signal weighting equation: weight_power[float]
restoring contrast: contrast [bool] (1 for denoising, 0 for brightening).
From https://github.com/BeataWereszczynska/k-space_wght_msk_for_MRI_denoising/blob/main/wght_msk_kspace.py
"""
# k-space weighting
if contrast:
kspace = kspace * np.power(abs(kspace), 3*weight_power)
# contrasting
a = kspace[abs(kspace) > abs(np.max(kspace))/6]
a = a / np.power(abs(a), weight_power)
kspace[abs(kspace) > abs(np.max(kspace))/6] = a
del a
else:
kspace = kspace * np.power(abs(kspace), weight_power)
del contrast, weight_power
# k-space masking
r = int(kspace.shape[0]/2)
mask = np.zeros(shape=kspace.shape)
cv2.circle(img=mask, center=(r,r), radius = r, color =(1,0,0), thickness=-1)
kspace = np.multiply(kspace, mask)
# return data
return kspace
def grad_mask_kspace(kspace, r):
"""
Graduate k-space masking for MRI image blurring.
Input:
k-space data: kspace [array],
radius for k-space masking in pixels: r [int].
From https://github.com/BeataWereszczynska/k-space_masking_for_MRI_denoising/blob/main/grad_mask_kspace.py
"""
# graduate k-space masking
mask_denoise = np.zeros(shape=kspace.shape)
mask = np.zeros(shape=kspace.shape)
for value in range(r, r+15, 2):
cv2.circle(img=mask, center=(int(kspace.shape[0]/2),int(kspace.shape[1]/2)),
radius = value, color =(1,0,0), thickness=-1)
mask_denoise = mask_denoise + mask
kspace = np.multiply(kspace, mask_denoise)
# return data
return kspace
def FFT_2D(kspace):
"""
Performs two-dimentional Fast Fourier Transform on k-space to reconstruct MRI image
(for Agilent FID data).
"""
# reconstruct MRI image
ft = np.fft.fft2(kspace) # 2D FFT
ft = np.fft.fftshift(ft) # fixing problem with corner being center of the image
ft = np.transpose(np.flip(ft, (1,0))) # matching geometry with VnmrJ-calculated image (still a bit shifted)
# return data
return ft
def kspace_det_edg(kspace, radius_min, radius_max, radius_step, threshold):
"""
k-space based details/edges detection in MRI images with optional k-space based denoising
and detail control.
Input:
k-space data: kspace[array],
denoising option: denoise [tuple (1=on 0=off [bool], weight_power [float], contrast [bool])
radii for k-space masking in pixels: radii [tuple (min [int], max [int], step [int])],
detail control option: detail_contr [tuple (1=on 0=off [bool], radius in pixels [int])]
threshold option: threshold [tuple (type [str], value [int])].
Threshold options are: ('auto', ), ('manual', threshold value), ('adaptive', C value)
"""
# masking k-space for creating an image of k-space-defined details/edges
masks = np.zeros(shape=kspace.shape)
mask = np.ones(shape=kspace.shape)
for value in range(radius_min, radius_max+1, radius_step):
cv2.circle(img=mask, center=(int(kspace.shape[0]/2),int(kspace.shape[1]/2)),
radius = value, color =(0,0,0), thickness=-1)
masks = masks + mask
kspace = np.multiply(kspace, masks)
# reconstructing the image of k-space-defined details/edges
ft2 = FFT_2D(kspace)
del kspace
# normalizing the image of k-space-defined details/edges (0-255)
ft2 = cv2.normalize(abs(ft2), None, 0, 255, cv2.NORM_MINMAX, dtype=cv2.CV_8U)
# k-space-defined details/edges to binary image
if threshold[0] == 'adaptive':
im_bw = cv2.adaptiveThreshold(ft2, 255, cv2.ADAPTIVE_THRESH_GAUSSIAN_C,
cv2.THRESH_BINARY, 31, threshold[1])
elif threshold[0] == 'manual':
im_bw = cv2.threshold(ft2, threshold[1], 255, cv2.THRESH_BINARY)[1]
elif threshold[0] == 'auto':
im_bw = cv2.threshold(ft2,0,255,cv2.THRESH_BINARY+cv2.THRESH_OTSU)[1]
else:
print("Threshold options are: ('auto', ), ('manual', threshold value), ('adaptive', C value).")
im_bw = np.zeros(ft2.shape)
del threshold
# denoise the binary image (remove small objects)
im_bw = morphology.remove_small_objects(im_bw != 0, min_size=3, in_place=True, connectivity=2)
# return data
return ft2, im_bw
def visualize(ft1, ft2, im_bw):
"""
Creates summarising figure.
"""
plt.rcParams['figure.dpi'] = 1200
plt.subplot(131)
plt.title('MRI image', fontdict = {'fontsize' : 8}), plt.axis('off')
plt.imshow(abs(ft1), cmap=plt.get_cmap('gray'))
plt.subplot(132)
plt.title('k-space-derived details/edges', fontdict = {'fontsize' : 8}), plt.axis('off')
plt.imshow(ft2, cmap=plt.get_cmap('gray'))
plt.subplot(133)
plt.title('binary image of details/edges', fontdict = {'fontsize' : 8}), plt.axis('off')
plt.imshow(im_bw, cmap=plt.get_cmap('gray'))
plt.tight_layout(pad=0, w_pad=0.2, h_pad=1.0)
plt.show()
def main():
# import parameters
path = 'sems_20190203_03.fid' # .fid folder location [str]
number_of_slices = 6 # total number of slices in the MRI experiment [int]
picked_slice = 4 # selected slice number [int]
# denoising parameters
denoise = 1 # 1=on 0=off [bool]
weight_power = 0.02 # exponent in the signal weighting equation [float]
contrast = 0 # restoring contrast: 1=on 0=off [bool]
# detail control parameters
detail_contr = 1 # 1=on 0=off [bool]
detail_r = 200 # radius for graduate masking in pixels [int]
# details/edges detection parameters
radius_min = 4 # smalest radius for k-space masking [int]
radius_max = 40 # largest radius for k-space masking [int]
radius_step = 1 # step betwen the two above values [int]
threshold = ('adaptive', -18) # threshold option [tuple (type, value)]
# threshold options: ('auto', ), ('manual', threshold value), ('adaptive', C value)
# import
kspace = import_kspace(path, number_of_slices, picked_slice)
del path, number_of_slices, picked_slice
# reconstructing the original MRI image
ft1 = FFT_2D(kspace)
# denoising
if denoise:
kspace = wght_msk_kspace(kspace, weight_power, contrast)
del weight_power, contrast, denoise
# detail control
if detail_contr:
kspace = grad_mask_kspace(kspace, detail_r)
del detail_r, detail_contr
# creating the image of details/edges and its binary version
ft2, im_bw = kspace_det_edg(kspace, radius_min, radius_max, radius_step, threshold)
del kspace
# visualization
visualize(ft1, ft2, im_bw)
# creating global variables to be available after the run completion
global MRI_complex_img
MRI_complex_img = ft1
global details_img
details_img = ft2
global binary_details
binary_details = im_bw
if __name__ == "__main__":
main()