├── GRAPPA_Tutorial.m ├── README.md └── fig ├── 01_SLphantom.png ├── 02_CoilChSensitivities.png ├── 03_SLphantomByCh.png ├── 04_kSpaceFullySampled.png ├── 05_kSpaceUnderSampled.png ├── 06_SLphantomAliased.png ├── 07_AutocalibrationLines.png ├── 08_kSpaceRestored.png ├── 09_ReconstructedImageByChannel.png └── 10_ReconstructedImage.png /GRAPPA_Tutorial.m: -------------------------------------------------------------------------------- 1 | % Code by Tetiana Dadakova 2 | % 04.21.2018 3 | % tetiana.d@gmail.com 4 | 5 | % This is a tutorial on MR image reconstruction using GRAPPA 6 | % Shepp-Logan phantom is first undersampled and then reconstructed using 7 | % GRAPPA 8 | 9 | close all 10 | 11 | %% 1. Load and display Shepp-Logan phantom 12 | 13 | phantom_temp = phantom('Modified Shepp-Logan',64); 14 | % Create complex phantom image 15 | phantom_shl = phantom_temp .* exp(1i * phantom_temp); 16 | 17 | % Display magnitude, phase, real, and imaginary parts of the phantom image 18 | figure; 19 | subplot(2,2,1); imshow(abs(phantom_shl),[]); title('abs') 20 | subplot(2,2,2); imshow(angle(phantom_shl),[]); title('angle') 21 | subplot(2,2,3); imshow(real(phantom_shl),[]); title('real') 22 | subplot(2,2,4); imshow(imag(phantom_shl),[]); title('imag') 23 | 24 | %% 2. Create phantom images for artificial coil channel sensitivities 25 | 26 | % Artificial coil channels sensitivities (6 channels): 27 | % linear gradients in 6 directions 28 | % adding more channels will improve the reconstruction 29 | ramp1 = double(repmat(1:2:128, [64, 1])); 30 | ramp2 = double(repmat(1:2:256, [128, 1])); 31 | ch_sensitivity_1 = ramp1/64; % left-right gradient 32 | ch_sensitivity_2 = imrotate(ch_sensitivity_1, 90); % bottom-up gradient 33 | ch_sensitivity_3 = imrotate(ch_sensitivity_1, 180); % right-left gradient 34 | ch_sensitivity_4 = imrotate(ch_sensitivity_1, 270); % up-bottom gradient 35 | temp = imrotate(ramp2/64, 45); 36 | ch_sensitivity_5 = temp(60:123,60:123); % diagonal left-right gradient 37 | ch_sensitivity_6 = imrotate(ch_sensitivity_5, 180); % diagonal right-left 38 | 39 | % Display channels sensitivities 40 | figure; 41 | subplot(2,3,1);imshow(ch_sensitivity_1,[]); title('Ch 1 sensitivity') 42 | subplot(2,3,2);imshow(ch_sensitivity_2,[]); title('Ch 2 sensitivity') 43 | subplot(2,3,3);imshow(ch_sensitivity_3,[]); title('Ch 3 sensitivity') 44 | subplot(2,3,4);imshow(ch_sensitivity_4,[]); title('Ch 4 sensitivity') 45 | subplot(2,3,5);imshow(ch_sensitivity_5,[]); title('Ch 5 sensitivity') 46 | subplot(2,3,6);imshow(ch_sensitivity_6,[]); title('Ch 6 sensitivity') 47 | 48 | % Phantom image accounted for coil channels sensitivities 49 | phantom_ch_1 = phantom_shl .* ch_sensitivity_1; 50 | phantom_ch_2 = phantom_shl .* ch_sensitivity_2; 51 | phantom_ch_3 = phantom_shl .* ch_sensitivity_3; 52 | phantom_ch_4 = phantom_shl .* ch_sensitivity_4; 53 | phantom_ch_5 = phantom_shl .* ch_sensitivity_5; 54 | phantom_ch_6 = phantom_shl .* ch_sensitivity_6; 55 | 56 | % Display magnitude 57 | figure; 58 | subplot(2,3,1);imshow(abs(phantom_ch_1),[]); 59 | subplot(2,3,2);imshow(abs(phantom_ch_2),[]); 60 | subplot(2,3,3);imshow(abs(phantom_ch_3),[]); 61 | subplot(2,3,4);imshow(abs(phantom_ch_4),[]); 62 | subplot(2,3,5);imshow(abs(phantom_ch_5),[]); 63 | subplot(2,3,6);imshow(abs(phantom_ch_6),[]); 64 | 65 | %% 3. Create undersampled k-space 66 | 67 | % 2D Fourier transform image to frequency space 68 | phantom_ch_1_k = fftshift(fft(ifftshift(phantom_ch_1,1),[],1),1); 69 | phantom_ch_1_k = fftshift(fft(ifftshift(phantom_ch_1_k,2),[],2),2); 70 | phantom_ch_2_k = fftshift(fft(ifftshift(phantom_ch_2,1),[],1),1); 71 | phantom_ch_2_k = fftshift(fft(ifftshift(phantom_ch_2_k,2),[],2),2); 72 | phantom_ch_3_k = fftshift(fft(ifftshift(phantom_ch_3,1),[],1),1); 73 | phantom_ch_3_k = fftshift(fft(ifftshift(phantom_ch_3_k,2),[],2),2); 74 | phantom_ch_4_k = fftshift(fft(ifftshift(phantom_ch_4,1),[],1),1); 75 | phantom_ch_4_k = fftshift(fft(ifftshift(phantom_ch_4_k,2),[],2),2); 76 | phantom_ch_5_k = fftshift(fft(ifftshift(phantom_ch_5,1),[],1),1); 77 | phantom_ch_5_k = fftshift(fft(ifftshift(phantom_ch_5_k,2),[],2),2); 78 | phantom_ch_6_k = fftshift(fft(ifftshift(phantom_ch_6,1),[],1),1); 79 | phantom_ch_6_k = fftshift(fft(ifftshift(phantom_ch_6_k,2),[],2),2); 80 | 81 | % Display fully-sampled k-space 82 | figure; 83 | subplot(2,3,1); imshow(abs(phantom_ch_1_k),[1 50]); title('Ch 1 k-space') 84 | subplot(2,3,2); imshow(abs(phantom_ch_2_k),[1 50]); title('Ch 2 k-space') 85 | subplot(2,3,3); imshow(abs(phantom_ch_3_k),[1 50]); title('Ch 3 k-space') 86 | subplot(2,3,4); imshow(abs(phantom_ch_4_k),[1 50]); title('Ch 4 k-space') 87 | subplot(2,3,5); imshow(abs(phantom_ch_5_k),[1 50]); title('Ch 5 k-space') 88 | subplot(2,3,6); imshow(abs(phantom_ch_6_k),[1 50]); title('Ch 6 k-space') 89 | 90 | % k-space undersampled twice (each second line is set to zeros) 91 | phantom_ch_1_k_u = zeros(size(phantom_ch_1_k)); 92 | phantom_ch_1_k_u(1:2:end,:) = phantom_ch_1_k(1:2:end,:); 93 | phantom_ch_2_k_u = zeros(size(phantom_ch_2_k)); 94 | phantom_ch_2_k_u(1:2:end,:) = phantom_ch_2_k(1:2:end,:); 95 | phantom_ch_3_k_u = zeros(size(phantom_ch_3_k)); 96 | phantom_ch_3_k_u(1:2:end,:) = phantom_ch_3_k(1:2:end,:); 97 | phantom_ch_4_k_u = zeros(size(phantom_ch_4_k)); 98 | phantom_ch_4_k_u(1:2:end,:) = phantom_ch_4_k(1:2:end,:); 99 | phantom_ch_5_k_u = zeros(size(phantom_ch_5_k)); 100 | phantom_ch_5_k_u(1:2:end,:) = phantom_ch_5_k(1:2:end,:); 101 | phantom_ch_6_k_u = zeros(size(phantom_ch_6_k)); 102 | phantom_ch_6_k_u(1:2:end,:) = phantom_ch_6_k(1:2:end,:); 103 | 104 | % Display undersampled k-space 105 | figure; 106 | subplot(2,3,1); imshow(abs(phantom_ch_1_k_u),[1 50]); title('Ch 1 undersampled k-space') 107 | subplot(2,3,2); imshow(abs(phantom_ch_2_k_u),[1 50]); title('Ch 2 undersampled k-space') 108 | subplot(2,3,3); imshow(abs(phantom_ch_3_k_u),[1 50]); title('Ch 3 undersampled k-space') 109 | subplot(2,3,4); imshow(abs(phantom_ch_4_k_u),[1 50]); title('Ch 4 undersampled k-space') 110 | subplot(2,3,5); imshow(abs(phantom_ch_5_k_u),[1 50]); title('Ch 5 undersampled k-space') 111 | subplot(2,3,6); imshow(abs(phantom_ch_6_k_u),[1 50]); title('Ch 6 undersampled k-space') 112 | 113 | % FT back to image space, to illustrate the ghosting artifact 114 | phantom_ch_1_im_u = fftshift(ifft(ifftshift(phantom_ch_1_k_u,1),[],1),1); 115 | phantom_ch_1_im_u = fftshift(ifft(ifftshift(phantom_ch_1_im_u,2),[],2),2); 116 | phantom_ch_2_im_u = fftshift(ifft(ifftshift(phantom_ch_2_k_u,1),[],1),1); 117 | phantom_ch_2_im_u = fftshift(ifft(ifftshift(phantom_ch_2_im_u,2),[],2),2); 118 | phantom_ch_3_im_u = fftshift(ifft(ifftshift(phantom_ch_3_k_u,1),[],1),1); 119 | phantom_ch_3_im_u = fftshift(ifft(ifftshift(phantom_ch_3_im_u,2),[],2),2); 120 | phantom_ch_4_im_u = fftshift(ifft(ifftshift(phantom_ch_4_k_u,1),[],1),1); 121 | phantom_ch_4_im_u = fftshift(ifft(ifftshift(phantom_ch_4_im_u,2),[],2),2); 122 | phantom_ch_5_im_u = fftshift(ifft(ifftshift(phantom_ch_5_k_u,1),[],1),1); 123 | phantom_ch_5_im_u = fftshift(ifft(ifftshift(phantom_ch_5_im_u,2),[],2),2); 124 | phantom_ch_6_im_u = fftshift(ifft(ifftshift(phantom_ch_6_k_u,1),[],1),1); 125 | phantom_ch_6_im_u = fftshift(ifft(ifftshift(phantom_ch_6_im_u,2),[],2),2); 126 | 127 | phantom_u_magn = abs(sqrt(phantom_ch_1_im_u.^2 + phantom_ch_2_im_u.^2 + ... 128 | phantom_ch_3_im_u.^2 + phantom_ch_4_im_u.^2 + ... 129 | phantom_ch_5_im_u.^2 + phantom_ch_6_im_u.^2)); 130 | figure; 131 | imshow(imresize(phantom_u_magn,5),[]); 132 | title({'Magnitude image of phantom', 'from twice undersampled k-space'}) 133 | 134 | %% 4. Choose autocalibration lines in kspace 135 | 136 | % Choose middle 16 lines in k-space for autocalibration 137 | phantom_ch_1_k_acl = zeros(size(phantom_ch_1_k)); 138 | phantom_ch_1_k_acl(25:40,:) = phantom_ch_1_k(25:40,:); 139 | phantom_ch_2_k_acl = zeros(size(phantom_ch_2_k)); 140 | phantom_ch_2_k_acl(25:40,:) = phantom_ch_2_k(25:40,:); 141 | phantom_ch_3_k_acl = zeros(size(phantom_ch_3_k)); 142 | phantom_ch_3_k_acl(25:40,:) = phantom_ch_3_k(25:40,:); 143 | phantom_ch_4_k_acl = zeros(size(phantom_ch_4_k)); 144 | phantom_ch_4_k_acl(25:40,:) = phantom_ch_4_k(25:40,:); 145 | phantom_ch_5_k_acl = zeros(size(phantom_ch_5_k)); 146 | phantom_ch_5_k_acl(25:40,:) = phantom_ch_5_k(25:40,:); 147 | phantom_ch_6_k_acl = zeros(size(phantom_ch_6_k)); 148 | phantom_ch_6_k_acl(25:40,:) = phantom_ch_6_k(25:40,:); 149 | 150 | % Display them 151 | figure; 152 | subplot(2,3,1); imshow(abs(phantom_ch_1_k_acl),[1 50]); title('Ch 1 autocalibration lines') 153 | subplot(2,3,2); imshow(abs(phantom_ch_2_k_acl),[1 50]); title('Ch 2 autocalibration lines') 154 | subplot(2,3,3); imshow(abs(phantom_ch_3_k_acl),[1 50]); title('Ch 3 autocalibration lines') 155 | subplot(2,3,4); imshow(abs(phantom_ch_4_k_acl),[1 50]); title('Ch 4 autocalibration lines') 156 | subplot(2,3,5); imshow(abs(phantom_ch_5_k_acl),[1 50]); title('Ch 5 autocalibration lines') 157 | subplot(2,3,6); imshow(abs(phantom_ch_6_k_acl),[1 50]); title('Ch 6 autocalibration lines') 158 | 159 | %% 5. Move kernel of size 3x3 to create the matrix of weights 160 | 161 | % Reshape each 3x3 patch in a vector and put all of them into a temp source 162 | % matrix for each channel of a coil 163 | kNo = 1; % kernel/patch number 164 | for ny = 25:(40-2) 165 | for nx = 1:(64-2) 166 | S_ch_1_temp(kNo,:) = ... 167 | reshape(phantom_ch_1_k_acl(ny:ny+2,nx:nx+2)',[1,9]); % ch1: each patch 3x3 is reshaped into vector and put into matrix one line after another 168 | S_ch_2_temp(kNo,:) = ... 169 | reshape(phantom_ch_2_k_acl(ny:ny+2,nx:nx+2)',[1,9]); 170 | S_ch_3_temp(kNo,:) = ... 171 | reshape(phantom_ch_3_k_acl(ny:ny+2,nx:nx+2)',[1,9]); 172 | S_ch_4_temp(kNo,:) = ... 173 | reshape(phantom_ch_4_k_acl(ny:ny+2,nx:nx+2)',[1,9]); 174 | S_ch_5_temp(kNo,:) = ... 175 | reshape(phantom_ch_5_k_acl(ny:ny+2,nx:nx+2)',[1,9]); 176 | S_ch_6_temp(kNo,:) = ... 177 | reshape(phantom_ch_6_k_acl(ny:ny+2,nx:nx+2)',[1,9]); 178 | kNo = kNo + 1; % to move through all patches 179 | end 180 | end 181 | % size(S_ch_1_temp) 182 | 183 | % Remove three middle ("unknown") values 184 | % The remaiming values form source matrix S, for each channel 185 | S_ch_1 = S_ch_1_temp(:,[1:3,7:9]); 186 | S_ch_2 = S_ch_2_temp(:,[1:3,7:9]); 187 | S_ch_3 = S_ch_3_temp(:,[1:3,7:9]); 188 | S_ch_4 = S_ch_4_temp(:,[1:3,7:9]); 189 | S_ch_5 = S_ch_5_temp(:,[1:3,7:9]); 190 | S_ch_6 = S_ch_6_temp(:,[1:3,7:9]); 191 | 192 | % Middle points form target vector T for each channel 193 | T_ch_1 = S_ch_1_temp(:,5); 194 | T_ch_2 = S_ch_2_temp(:,5); 195 | T_ch_3 = S_ch_3_temp(:,5); 196 | T_ch_4 = S_ch_4_temp(:,5); 197 | T_ch_5 = S_ch_5_temp(:,5); 198 | T_ch_6 = S_ch_6_temp(:,5); 199 | 200 | % Stack all the channels together in the same matrix one after another 201 | S = [S_ch_1 S_ch_2 S_ch_3 S_ch_4 S_ch_5 S_ch_6]; 202 | T = [T_ch_1 T_ch_2 T_ch_3 T_ch_4 T_ch_5 T_ch_6]; 203 | 204 | % Invert S to find weights 205 | W = pinv(S) * T; 206 | 207 | %% 6. Forward problem to find missing lines 208 | 209 | % Construct source matric from the undersampled image 210 | kNo = 1; % kernel/patch number 211 | for ny = 1:2:(64-2) 212 | for nx = 1:(64-2) 213 | S_ch_1_new_temp(kNo,:) = ... 214 | reshape(phantom_ch_1_k_u(ny:ny+2,nx:nx+2)',[1,9]); 215 | S_ch_2_new_temp(kNo,:) = ... 216 | reshape(phantom_ch_2_k_u(ny:ny+2,nx:nx+2)',[1,9]); 217 | S_ch_3_new_temp(kNo,:) = ... 218 | reshape(phantom_ch_3_k_u(ny:ny+2,nx:nx+2)',[1,9]); 219 | S_ch_4_new_temp(kNo,:) = ... 220 | reshape(phantom_ch_4_k_u(ny:ny+2,nx:nx+2)',[1,9]); 221 | S_ch_5_new_temp(kNo,:) = ... 222 | reshape(phantom_ch_5_k_u(ny:ny+2,nx:nx+2)',[1,9]); 223 | S_ch_6_new_temp(kNo,:) = ... 224 | reshape(phantom_ch_6_k_u(ny:ny+2,nx:nx+2)',[1,9]); 225 | kNo = kNo + 1; 226 | end 227 | end 228 | % size(S_ch_1_new_temp) 229 | 230 | % Remove three middle ("unknown") values 231 | % The remaiming values form matrix S, for each channel 232 | S_ch_1_new = S_ch_1_new_temp(:,[1:3,7:9]); 233 | S_ch_2_new = S_ch_2_new_temp(:,[1:3,7:9]); 234 | S_ch_3_new = S_ch_3_new_temp(:,[1:3,7:9]); 235 | S_ch_4_new = S_ch_4_new_temp(:,[1:3,7:9]); 236 | S_ch_5_new = S_ch_5_new_temp(:,[1:3,7:9]); 237 | S_ch_6_new = S_ch_6_new_temp(:,[1:3,7:9]); 238 | 239 | % Stack all the channels together in the same matrix one after another 240 | S_new = [S_ch_1_new S_ch_2_new S_ch_3_new S_ch_4_new S_ch_5_new S_ch_6_new]; 241 | 242 | % T_unknown = S_undersampled * W 243 | T_new = S_new * W; 244 | 245 | %% 7. Filling in the missing lines into undersampled image 246 | 247 | T_ch_1_new_M = reshape(T_new(:,1),[62,31]); % reshape vector into matrix 248 | T_ch_2_new_M = reshape(T_new(:,2),[62,31]); 249 | T_ch_3_new_M = reshape(T_new(:,3),[62,31]); 250 | T_ch_4_new_M = reshape(T_new(:,4),[62,31]); 251 | T_ch_5_new_M = reshape(T_new(:,5),[62,31]); 252 | T_ch_6_new_M = reshape(T_new(:,6),[62,31]); 253 | 254 | T_ch_1_new_M = T_ch_1_new_M'; 255 | T_ch_2_new_M = T_ch_2_new_M'; 256 | T_ch_3_new_M = T_ch_3_new_M'; 257 | T_ch_4_new_M = T_ch_4_new_M'; 258 | T_ch_5_new_M = T_ch_5_new_M'; 259 | T_ch_6_new_M = T_ch_6_new_M'; 260 | 261 | % Fill in approximated lines 262 | P1_f_u_new = phantom_ch_1_k_u; 263 | P2_f_u_new = phantom_ch_2_k_u; 264 | P3_f_u_new = phantom_ch_3_k_u; 265 | P4_f_u_new = phantom_ch_4_k_u; 266 | P5_f_u_new = phantom_ch_5_k_u; 267 | P6_f_u_new = phantom_ch_6_k_u; 268 | P1_f_u_new(2:2:end-1,2:end-1) = T_ch_1_new_M; 269 | P2_f_u_new(2:2:end-1,2:end-1) = T_ch_2_new_M; 270 | P3_f_u_new(2:2:end-1,2:end-1) = T_ch_3_new_M; 271 | P4_f_u_new(2:2:end-1,2:end-1) = T_ch_4_new_M; 272 | P5_f_u_new(2:2:end-1,2:end-1) = T_ch_5_new_M; 273 | P6_f_u_new(2:2:end-1,2:end-1) = T_ch_6_new_M; 274 | 275 | % Display undersampled k-space 276 | figure; 277 | subplot(2,3,1); imshow(abs(phantom_ch_1_k_u),[1 50]); title('Ch 1 undersampled k-space') 278 | subplot(2,3,2); imshow(abs(phantom_ch_2_k_u),[1 50]); title('Ch 2 undersampled k-space') 279 | subplot(2,3,3); imshow(abs(phantom_ch_3_k_u),[1 50]); title('Ch 3 undersampled k-space') 280 | subplot(2,3,4); imshow(abs(phantom_ch_4_k_u),[1 50]); title('Ch 4 undersampled k-space') 281 | subplot(2,3,5); imshow(abs(phantom_ch_5_k_u),[1 50]); title('Ch 5 undersampled k-space') 282 | subplot(2,3,6); imshow(abs(phantom_ch_6_k_u),[1 50]); title('Ch 6 undersampled k-space') 283 | 284 | % Display k-space after the approximated lines were filled in 285 | figure; 286 | subplot(2,3,1); imshow(abs(P1_f_u_new),[1 50]); title('Ch 1 restored k-space') 287 | subplot(2,3,2); imshow(abs(P2_f_u_new),[1 50]); title('Ch 2 restored k-space') 288 | subplot(2,3,3); imshow(abs(P3_f_u_new),[1 50]); title('Ch 3 restored k-space') 289 | subplot(2,3,4); imshow(abs(P4_f_u_new),[1 50]); title('Ch 4 restored k-space') 290 | subplot(2,3,5); imshow(abs(P5_f_u_new),[1 50]); title('Ch 5 restored k-space') 291 | subplot(2,3,6); imshow(abs(P6_f_u_new),[1 50]); title('Ch 6 restored k-space') 292 | 293 | % Fourier transform restored k-space for each channel 294 | Im_Recon_ch_1 = fftshift(ifft(ifftshift(P1_f_u_new,1),[],1),1); 295 | Im_Recon_ch_1 = fftshift(ifft(ifftshift(Im_Recon_ch_1,2),[],2),2); 296 | Im_Recon_ch_2 = fftshift(ifft(ifftshift(P2_f_u_new,1),[],1),1); 297 | Im_Recon_ch_2 = fftshift(ifft(ifftshift(Im_Recon_ch_2,2),[],2),2); 298 | Im_Recon_ch_3 = fftshift(ifft(ifftshift(P3_f_u_new,1),[],1),1); 299 | Im_Recon_ch_3 = fftshift(ifft(ifftshift(Im_Recon_ch_3,2),[],2),2); 300 | Im_Recon_ch_4 = fftshift(ifft(ifftshift(P4_f_u_new,1),[],1),1); 301 | Im_Recon_ch_4 = fftshift(ifft(ifftshift(Im_Recon_ch_4,2),[],2),2); 302 | Im_Recon_ch_5 = fftshift(ifft(ifftshift(P5_f_u_new,1),[],1),1); 303 | Im_Recon_ch_5 = fftshift(ifft(ifftshift(Im_Recon_ch_5,2),[],2),2); 304 | Im_Recon_ch_6 = fftshift(ifft(ifftshift(P6_f_u_new,1),[],1),1); 305 | Im_Recon_ch_6 = fftshift(ifft(ifftshift(Im_Recon_ch_6,2),[],2),2); 306 | 307 | % Display restores images for each channel 308 | figure; 309 | subplot(2,3,1); imshow(abs(Im_Recon_ch_1),[]); title('Reconstructed image, ch 1') 310 | subplot(2,3,2); imshow(abs(Im_Recon_ch_2),[]); title('Reconstructed image, ch 2') 311 | subplot(2,3,3); imshow(abs(Im_Recon_ch_3),[]); title('Reconstructed image, ch 3') 312 | subplot(2,3,4); imshow(abs(Im_Recon_ch_4),[]); title('Reconstructed image, ch 4') 313 | subplot(2,3,5); imshow(abs(Im_Recon_ch_5),[]); title('Reconstructed image, ch 5') 314 | subplot(2,3,6); imshow(abs(Im_Recon_ch_6),[]); title('Reconstructed image, ch 6') 315 | 316 | % Combine coil channels 317 | Im_Recon = sqrt(abs(Im_Recon_ch_1).^2 + abs(Im_Recon_ch_2).^2 + ... 318 | abs(Im_Recon_ch_3).^2 + abs(Im_Recon_ch_4).^2 + ... 319 | abs(Im_Recon_ch_5).^2 + abs(Im_Recon_ch_6).^2); 320 | 321 | % Display image reconstructed uisng GRAPPA 322 | figure; imshow(imresize(Im_Recon,5),[]); title('Image reconstructed using GRAPPA') 323 | 324 | 325 | 326 | 327 | -------------------------------------------------------------------------------- /README.md: -------------------------------------------------------------------------------- 1 | # Tutorial-MRI-Reconstruction-Using-GRAPPA 2 | I wrote this code when trying to understand how GRAPPA reconstruction works. I was trying to keep it very simple and readable. Thought, it might be useful for someone else. Let me know if you have questions, comments, or suggestions: tetiana.d@gmail.com 3 | 4 | 5 | Code consists of the following parts: 6 | 7 | 1. Load Shepp-Logan phantom. 8 | ![alt text](fig/01_SLphantom.png) 9 | 10 | 2. Create six artificial sensitivities of the channels of the coil. This is done by creating linear intensity gradients in 6 directions. Adding more channels will improve the reconstruction. 11 | ![alt text](fig/02_CoilChSensitivities.png) 12 | ![alt text](fig/03_SLphantomByCh.png) 13 | 14 | 3. Fourier transform each channel image to k-space 15 | ![alt text](fig/04_kSpaceFullySampled.png) 16 | and undersample it twice. 17 | ![alt text](fig/05_kSpaceUnderSampled.png) 18 | Reconstruct the Shepp-Logan phantom image from the under sampled k-space to show the artifact. 19 | ![alt text](fig/06_SLphantomAliased.png) 20 | 21 | 4. Choose auto calibration lines. 22 | ![alt text](fig/07_AutocalibrationLines.png) 23 | 24 | 5. Source (S) pixels are related to target (T) pixels by weighs (W): S * W = T. To find weights: W = inv(S) * T. (See schematic figure here: http://mriquestions.com/grappaarc.html) 25 | 26 | 6. Forward problem to find missing lines. T_unknown = S_undersampled * W 27 | 28 | 7. Fill in the missing lines into undersampled image. 29 | ![alt text](fig/08_kSpaceRestored.png) 30 | Reconstruct image for each channel. 31 | ![alt text](fig/09_ReconstructedImageByChannel.png) 32 | And finally reconstruct the image. 33 | ![alt text](fig/10_ReconstructedImage.png) 34 | 35 | 36 | 37 | 38 | -------------------------------------------------------------------------------- /fig/01_SLphantom.png: -------------------------------------------------------------------------------- https://raw.githubusercontent.com/tetianadadakova/Tutorial-MRI-Reconstruction-Using-GRAPPA/9976634988e79ec9baa718c3f3cfbab774596312/fig/01_SLphantom.png -------------------------------------------------------------------------------- /fig/02_CoilChSensitivities.png: -------------------------------------------------------------------------------- https://raw.githubusercontent.com/tetianadadakova/Tutorial-MRI-Reconstruction-Using-GRAPPA/9976634988e79ec9baa718c3f3cfbab774596312/fig/02_CoilChSensitivities.png -------------------------------------------------------------------------------- /fig/03_SLphantomByCh.png: -------------------------------------------------------------------------------- https://raw.githubusercontent.com/tetianadadakova/Tutorial-MRI-Reconstruction-Using-GRAPPA/9976634988e79ec9baa718c3f3cfbab774596312/fig/03_SLphantomByCh.png -------------------------------------------------------------------------------- /fig/04_kSpaceFullySampled.png: -------------------------------------------------------------------------------- https://raw.githubusercontent.com/tetianadadakova/Tutorial-MRI-Reconstruction-Using-GRAPPA/9976634988e79ec9baa718c3f3cfbab774596312/fig/04_kSpaceFullySampled.png -------------------------------------------------------------------------------- /fig/05_kSpaceUnderSampled.png: -------------------------------------------------------------------------------- https://raw.githubusercontent.com/tetianadadakova/Tutorial-MRI-Reconstruction-Using-GRAPPA/9976634988e79ec9baa718c3f3cfbab774596312/fig/05_kSpaceUnderSampled.png -------------------------------------------------------------------------------- /fig/06_SLphantomAliased.png: -------------------------------------------------------------------------------- https://raw.githubusercontent.com/tetianadadakova/Tutorial-MRI-Reconstruction-Using-GRAPPA/9976634988e79ec9baa718c3f3cfbab774596312/fig/06_SLphantomAliased.png -------------------------------------------------------------------------------- /fig/07_AutocalibrationLines.png: -------------------------------------------------------------------------------- https://raw.githubusercontent.com/tetianadadakova/Tutorial-MRI-Reconstruction-Using-GRAPPA/9976634988e79ec9baa718c3f3cfbab774596312/fig/07_AutocalibrationLines.png -------------------------------------------------------------------------------- /fig/08_kSpaceRestored.png: -------------------------------------------------------------------------------- https://raw.githubusercontent.com/tetianadadakova/Tutorial-MRI-Reconstruction-Using-GRAPPA/9976634988e79ec9baa718c3f3cfbab774596312/fig/08_kSpaceRestored.png -------------------------------------------------------------------------------- /fig/09_ReconstructedImageByChannel.png: -------------------------------------------------------------------------------- https://raw.githubusercontent.com/tetianadadakova/Tutorial-MRI-Reconstruction-Using-GRAPPA/9976634988e79ec9baa718c3f3cfbab774596312/fig/09_ReconstructedImageByChannel.png -------------------------------------------------------------------------------- /fig/10_ReconstructedImage.png: -------------------------------------------------------------------------------- https://raw.githubusercontent.com/tetianadadakova/Tutorial-MRI-Reconstruction-Using-GRAPPA/9976634988e79ec9baa718c3f3cfbab774596312/fig/10_ReconstructedImage.png --------------------------------------------------------------------------------