-
Notifications
You must be signed in to change notification settings - Fork 0
/
fanbeamtomo.m
309 lines (255 loc) · 8.97 KB
/
fanbeamtomo.m
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
function [A,b,x,theta,p,R,d] = fanbeamtomo(N,theta,p,R,d,isDisp, absCoeff)
%FANBEAMTOMO Creates a 2D tomography test problem using fan beams
%
% [A,b,x,theta,p,R,d] = fanbeamtomo(N)
% [A,b,x,theta,p,R,d] = fanbeamtomo(N,theta)
% [A,b,x,theta,p,R,d] = fanbeamtomo(N,theta,p)
% [A,b,x,theta,p,R,d] = fanbeamtomo(N,theta,p,R)
% [A,b,x,theta,p,R,d] = fanbeamtomo(N,theta,p,R,d)
% [A,b,x,theta,p,R,d] = fanbeamtomo(N,theta,p,R,d,isDisp)
%
% This function creates a 2D tomography test problem with an N-times-N
% domain, using p rays in fan-formation for each angle in the vector theta.
%
% Input:
% N Scalar denoting the number of discretization intervals in
% each dimesion, such that the domain consists of N^2 cells.
% theta Vector containing the angles in degrees. Default: theta =
% 0:1:359.
% p Number of rays for each angle. Default: p =
% round(sqrt(2)*N).
% R The distance from the source to the center of the domain
% is R*N. Default: R = 2.
% d Scalar that determines the span of the rays. The default
% value is defined such that from (0,R*N) the first ray hits
% the point (-N/2,N/2) and the last ray hits (N/2,N/2).
% isDisp If isDisp is non-zero it specifies the time in seconds
% to pause in the display of the rays. If zero (the default),
% no display is shown.
%
% Output:
% A Coefficient matrix with N^2 columns and nA*p rows,
% where nA is the number of angles, i.e., length(theta).
% b Vector containing the rhs of the test problem.
% x Vector containing the exact solution, with elements
% between 0 and 1.
% theta Vector containing the used angles in degrees.
% p The number of used rays for each angle.
% R The radius in side lengths.
% d The span of the rays.
%
% See also: paralleltomo, seismictomo.
% Jakob Sauer Jørgensen, Maria Saxild-Hansen and Per Christian Hansen,
% October 1, 2014, DTU Compute.
% Reference: A. C. Kak and M. Slaney, Principles of Computerized
% Tomographic Imaging, SIAM, Philadelphia, 2001.
% Default illustration:
if nargin < 6 || isempty(isDisp)
isDisp = 0;
end
% Default value of R.
if nargin < 4 || isempty(R)
R = 2;
end
% Default value of d.
if nargin < 5 || isempty(d)
% Determine angular span. The first and last rays touch points (-N/2,N/2)
% and (N/2,N/2), respectively.
d = 2*atand(1/(2*R-1));
end
% Default value of the number of rays p.
if nargin < 3 || isempty(p)
p = round(sqrt(2)*N);
end
% Default value of the angles theta.
if nargin < 2 || isempty(theta)
theta = 0:359;
end
% Input check. The source must lie outside the domain.
if R < sqrt(2)/2
error('R must be greater than half squareroot 2')
end
R = R*N;
% The width of the angle of the source.
if d < 0 || d > 180
error('The angle of the source must be in the interval [0 180]')
end
% Anonymous function rotation matrix
Omega_x = @(omega_par) [cosd(omega_par) -sind(omega_par)];
Omega_y = @(omega_par) [sind(omega_par) cosd(omega_par)];
% nA denotes the number of angles.
nA = length(theta);
% The starting values both the x and the y coordinates.
x0 = 0;
y0 = R;
xy0 = [x0;y0];
omega = linspace(-d/2,d/2,p);
% The intersection lines.
x = (-N/2:N/2)';
y = x;
% Initialize vectors that contains the row numbers, the column numbers and
% the values for creating the matrix A effiecently.
rows = zeros(2*N*nA*p,1);
cols = rows;
vals = rows;
idxend = 0;
% Prepare for illustration
if isDisp
AA = rand(N);
figure
end
% Loop over the chosen angles of the source.
for i = 1:nA
% The starting points for the current angle theta.
x0theta = Omega_x(theta(i))*xy0;
y0theta = Omega_y(theta(i))*xy0;
% The starting (center) direction vector (opposite xytheta) and
% normalized.
xytheta = [x0theta; y0theta];
abtheta = -xytheta/R;
% Illustration of the domain
if isDisp % illustration of source
clf
pause(isDisp)
imagesc((-N/2+.5):(N/2-0.5),(-N/2+.5):(N/2-0.5),AA), colormap gray,
axis xy
hold on
axis(1.1*R*[-1,1,-1,1])
axis equal
plot(x0theta,y0theta,'o','color',[60 179 113]/255,...
'linewidth',1.5,'markersize',10)
end
% Loop over the rays.
for j = 1:p
% The direction vector for the current ray.
a = Omega_x(omega(j))*abtheta;
b = Omega_y(omega(j))*abtheta;
% Illustration of rays
if isDisp
plot([x0theta,x0theta+1.7*R*a],[y0theta,y0theta+1.7*R*b],'-',...
'color',[220 0 0]/255,'linewidth',1.5)
axis(R*[-1,1,-1,1])
end
% Use the parametrisation of line to get the y-coordinates of
% intersections with x = k, i.e. x constant.
tx = (x - x0theta)/a;
yx = b*tx + y0theta;
% Use the parametrisation of line to get the x-coordinates of
% intersections with y = k, i.e. y constant.
ty = (y - y0theta)/b;
xy = a*ty + x0theta;
% Collect the intersection times and coordinates.
t = [tx; ty];
xxy = [x; xy];
yxy = [yx; y];
% Sort the coordinates according to intersection time.
[~,I] = sort(t);
xxy = xxy(I);
yxy = yxy(I);
% Skip the points outside the box.
I = (xxy >= -N/2 & xxy <= N/2 & yxy >= -N/2 & yxy <= N/2);
xxy = xxy(I);
yxy = yxy(I);
% Skip double points.
I = (abs(diff(xxy)) <= 1e-10 & abs(diff(yxy)) <= 1e-10);
xxy(I) = [];
yxy(I) = [];
% Illustration of the rays
if isDisp
set(gca,'Xticklabel',{})
set(gca,'Yticklabel',{})
pause(isDisp)
end
% Calculate the length within cell and determines the number of
% cells which is hit.
d = sqrt(diff(xxy).^2 + diff(yxy).^2);
numvals = numel(d);
% Store the values inside the box.
if numvals > 0
% Calculates the midpoints of the line within the cells.
xm = 0.5*(xxy(1:end-1)+xxy(2:end)) + N/2;
ym = 0.5*(yxy(1:end-1)+yxy(2:end)) + N/2;
% Translate the midpoint coordinates to index.
col = floor(xm)*N + (N - floor(ym));
% Create the indices to store the values to vector for
% later creation of A matrix.
idxstart = idxend + 1;
idxend = idxstart + numvals - 1;
idx = idxstart:idxend;
% Store row numbers, column numbers and values.
rows(idx) = (i-1)*p + j;
cols(idx) = col;
vals(idx) = d;
end
end
end
% Truncate excess zeros.
rows = rows(1:idxend);
cols = cols(1:idxend);
vals = vals(1:idxend);
% Create sparse matrix A from the stored values.
A = sparse(rows,cols,vals,p*nA,N^2);
if nargout > 1
% Create phantom head as a reshaped vector.
x = absCoeff;
x = x(:);
% Create rhs.
b = A*x;
end
if nargout > 5
R = R/N;
end
if nargout > 6
d = 2*omega(end); % Restore value of d.
end;
function X = myphantom(N)
%MYPHANTOM creates the modified Shepp-Logan phantom
% X = myphantom(N)
%
% This function create the modifed Shepp-Logan phantom with the
% discretization N x N, and returns it as a vector.
%
% Input:
% N Scalar denoting the nubmer of discretization intervals in each
% dimesion, such that the phantom head consists of N^2 cells.
%
% Output:
% X The modified phantom head reshaped as a vector
% This head phantom is the same as the Shepp-Logan except the intensities
% are changed to yield higher contrast in the image.
%
% Peter Toft, "The Radon Transform - Theory and Implementation", PhD
% thesis, DTU Informatics, Technical University of Denmark, June 1996.
% A a b x0 y0 phi
% ---------------------------------
e = [ 1 .69 .92 0 0 0
-.8 .6624 .8740 0 -.0184 0
-.2 .1100 .3100 .22 0 -18
-.2 .1600 .4100 -.22 0 18
.1 .2100 .2500 0 .35 0
.1 .0460 .0460 0 .1 0
.1 .0460 .0460 0 -.1 0
.1 .0460 .0230 -.08 -.605 0
.1 .0230 .0230 0 -.606 0
.1 .0230 .0460 .06 -.605 0 ];
xn = ((0:N-1)-(N-1)/2)/((N-1)/2);
Xn = repmat(xn,N,1);
Yn = rot90(Xn);
X = zeros(N);
% For each ellipse to be added
for i = 1:size(e,1)
a2 = e(i,2)^2;
b2 = e(i,3)^2;
x0 = e(i,4);
y0 = e(i,5);
phi = e(i,6)*pi/180;
A = e(i,1);
x = Xn-x0;
y = Yn-y0;
index = find(((x.*cos(phi) + y.*sin(phi)).^2)./a2 + ...
((y.*cos(phi) - x.*sin(phi))).^2./b2 <= 1);
% Add the amplitude of the ellipse
X(index) = X(index) + A;
end
% Return as vector and ensure nonnegative elements.
X = X(:); X(X<0) = 0;