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
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
|
/*
-----------------------------------------------------------------------
Copyright: 2010-2018, imec Vision Lab, University of Antwerp
2014-2018, CWI, Amsterdam
Contact: astra@astra-toolbox.com
Website: http://www.astra-toolbox.com/
This file is part of the ASTRA Toolbox.
The ASTRA Toolbox is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
The ASTRA Toolbox is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with the ASTRA Toolbox. If not, see <http://www.gnu.org/licenses/>.
-----------------------------------------------------------------------
*/
#include "astra/cuda/2d/util.h"
#include "astra/cuda/2d/arith.h"
#ifdef STANDALONE
#include "testutil.h"
#endif
#include <cstdio>
#include <cassert>
#include <iostream>
#include <list>
#include <cmath>
typedef texture<float, 2, cudaReadModeElementType> texture2D;
static texture2D gT_volumeTexture;
namespace astraCUDA {
static const unsigned g_MaxAngles = 2560;
__constant__ float gC_angle[g_MaxAngles];
__constant__ float gC_angle_offset[g_MaxAngles];
__constant__ float gC_angle_detsize[g_MaxAngles];
// optimization parameters
static const unsigned int g_anglesPerBlock = 16;
static const unsigned int g_detBlockSize = 32;
static const unsigned int g_blockSlices = 64;
// fixed point scaling factor
#define fPREC_FACTOR 16.0f
#define iPREC_FACTOR 16
static bool bindVolumeDataTexture(float* data, cudaArray*& dataArray, unsigned int pitch, unsigned int width, unsigned int height)
{
cudaChannelFormatDesc channelDesc = cudaCreateChannelDesc<float>();
dataArray = 0;
cudaMallocArray(&dataArray, &channelDesc, width, height);
cudaMemcpy2DToArray(dataArray, 0, 0, data, pitch*sizeof(float), width*sizeof(float), height, cudaMemcpyDeviceToDevice);
gT_volumeTexture.addressMode[0] = cudaAddressModeBorder;
gT_volumeTexture.addressMode[1] = cudaAddressModeBorder;
gT_volumeTexture.filterMode = cudaFilterModeLinear;
gT_volumeTexture.normalized = false;
// TODO: For very small sizes (roughly <=512x128) with few angles (<=180)
// not using an array is more efficient.
// cudaBindTexture2D(0, gT_volumeTexture, (const void*)data, channelDesc, width, height, sizeof(float)*pitch);
cudaBindTextureToArray(gT_volumeTexture, dataArray, channelDesc);
// TODO: error value?
return true;
}
// projection for angles that are roughly horizontal
// (detector roughly vertical)
__global__ void FPhorizontal_simple(float* D_projData, unsigned int projPitch, unsigned int startSlice, unsigned int startAngle, unsigned int endAngle, const SDimensions dims, float outputScale)
{
const int relDet = threadIdx.x;
const int relAngle = threadIdx.y;
int angle = startAngle + blockIdx.x * g_anglesPerBlock + relAngle;
if (angle >= endAngle)
return;
const float theta = gC_angle[angle];
const float cos_theta = __cosf(theta);
const float sin_theta = __sinf(theta);
// compute start detector for this block/angle:
const int detRegion = blockIdx.y;
const int detector = detRegion * g_detBlockSize + relDet;
// Now project the part of the ray to angle,detector through
// slices startSlice to startSlice+g_blockSlices-1
if (detector < 0 || detector >= dims.iProjDets)
return;
const float fDetStep = -gC_angle_detsize[angle] / sin_theta;
float fSliceStep = cos_theta / sin_theta;
float fDistCorr;
if (sin_theta > 0.0f)
fDistCorr = outputScale / sin_theta;
else
fDistCorr = -outputScale / sin_theta;
float fVal = 0.0f;
// project detector on slice
float fP = (detector - 0.5f*dims.iProjDets + 0.5f - gC_angle_offset[angle]) * fDetStep + (startSlice - 0.5f*dims.iVolWidth + 0.5f) * fSliceStep + 0.5f*dims.iVolHeight - 0.5f + 0.5f;
float fS = startSlice + 0.5f;
int endSlice = startSlice + g_blockSlices;
if (endSlice > dims.iVolWidth)
endSlice = dims.iVolWidth;
if (dims.iRaysPerDet > 1) {
fP += (-0.5f*dims.iRaysPerDet + 0.5f)/dims.iRaysPerDet * fDetStep;
const float fSubDetStep = fDetStep / dims.iRaysPerDet;
fDistCorr /= dims.iRaysPerDet;
fSliceStep -= dims.iRaysPerDet * fSubDetStep;
for (int slice = startSlice; slice < endSlice; ++slice)
{
for (int iSubT = 0; iSubT < dims.iRaysPerDet; ++iSubT) {
fVal += tex2D(gT_volumeTexture, fS, fP);
fP += fSubDetStep;
}
fP += fSliceStep;
fS += 1.0f;
}
} else {
for (int slice = startSlice; slice < endSlice; ++slice)
{
fVal += tex2D(gT_volumeTexture, fS, fP);
fP += fSliceStep;
fS += 1.0f;
}
}
D_projData[angle*projPitch+detector] += fVal * fDistCorr;
}
// projection for angles that are roughly vertical
// (detector roughly horizontal)
__global__ void FPvertical_simple(float* D_projData, unsigned int projPitch, unsigned int startSlice, unsigned int startAngle, unsigned int endAngle, const SDimensions dims, float outputScale)
{
const int relDet = threadIdx.x;
const int relAngle = threadIdx.y;
int angle = startAngle + blockIdx.x * g_anglesPerBlock + relAngle;
if (angle >= endAngle)
return;
const float theta = gC_angle[angle];
const float cos_theta = __cosf(theta);
const float sin_theta = __sinf(theta);
// compute start detector for this block/angle:
const int detRegion = blockIdx.y;
const int detector = detRegion * g_detBlockSize + relDet;
// Now project the part of the ray to angle,detector through
// slices startSlice to startSlice+g_blockSlices-1
if (detector < 0 || detector >= dims.iProjDets)
return;
const float fDetStep = gC_angle_detsize[angle] / cos_theta;
float fSliceStep = sin_theta / cos_theta;
float fDistCorr;
if (cos_theta < 0.0f)
fDistCorr = -outputScale / cos_theta;
else
fDistCorr = outputScale / cos_theta;
float fVal = 0.0f;
float fP = (detector - 0.5f*dims.iProjDets + 0.5f - gC_angle_offset[angle]) * fDetStep + (startSlice - 0.5f*dims.iVolHeight + 0.5f) * fSliceStep + 0.5f*dims.iVolWidth - 0.5f + 0.5f;
float fS = startSlice+0.5f;
int endSlice = startSlice + g_blockSlices;
if (endSlice > dims.iVolHeight)
endSlice = dims.iVolHeight;
if (dims.iRaysPerDet > 1) {
fP += (-0.5f*dims.iRaysPerDet + 0.5f)/dims.iRaysPerDet * fDetStep;
const float fSubDetStep = fDetStep / dims.iRaysPerDet;
fDistCorr /= dims.iRaysPerDet;
fSliceStep -= dims.iRaysPerDet * fSubDetStep;
for (int slice = startSlice; slice < endSlice; ++slice)
{
for (int iSubT = 0; iSubT < dims.iRaysPerDet; ++iSubT) {
fVal += tex2D(gT_volumeTexture, fP, fS);
fP += fSubDetStep;
}
fP += fSliceStep;
fS += 1.0f;
}
} else {
for (int slice = startSlice; slice < endSlice; ++slice)
{
fVal += tex2D(gT_volumeTexture, fP, fS);
fP += fSliceStep;
fS += 1.0f;
}
}
D_projData[angle*projPitch+detector] += fVal * fDistCorr;
}
// Coordinates of center of detector pixel number t:
// x = (t - 0.5*nDets + 0.5 - fOffset) * fSize * cos(fAngle)
// y = - (t - 0.5*nDets + 0.5 - fOffset) * fSize * sin(fAngle)
static void convertAndUploadAngles(const SParProjection *projs, unsigned int nth, unsigned int ndets)
{
float *angles = new float[nth];
float *offsets = new float[nth];
float *detsizes = new float[nth];
for (int i = 0; i < nth; ++i)
getParParameters(projs[i], ndets, angles[i], detsizes[i], offsets[i]);
cudaMemcpyToSymbol(gC_angle, angles, nth*sizeof(float), 0, cudaMemcpyHostToDevice);
cudaMemcpyToSymbol(gC_angle_offset, offsets, nth*sizeof(float), 0, cudaMemcpyHostToDevice);
cudaMemcpyToSymbol(gC_angle_detsize, detsizes, nth*sizeof(float), 0, cudaMemcpyHostToDevice);
}
bool FP_simple_internal(float* D_volumeData, unsigned int volumePitch,
float* D_projData, unsigned int projPitch,
const SDimensions& dims, const SParProjection* angles,
float outputScale)
{
assert(dims.iProjAngles <= g_MaxAngles);
assert(angles);
cudaArray* D_dataArray;
bindVolumeDataTexture(D_volumeData, D_dataArray, volumePitch, dims.iVolWidth, dims.iVolHeight);
convertAndUploadAngles(angles, dims.iProjAngles, dims.iProjDets);
dim3 dimBlock(g_detBlockSize, g_anglesPerBlock); // detector block size, angles
std::list<cudaStream_t> streams;
// Run over all angles, grouping them into groups of the same
// orientation (roughly horizontal vs. roughly vertical).
// Start a stream of grids for each such group.
// TODO: Check if it's worth it to store this info instead
// of recomputing it every FP.
unsigned int blockStart = 0;
unsigned int blockEnd = 0;
bool blockVertical = false;
for (unsigned int a = 0; a <= dims.iProjAngles; ++a) {
bool vertical = false;
// TODO: Having <= instead of < below causes a 5% speedup.
// Maybe we should detect corner cases and put them in the optimal
// group of angles.
if (a != dims.iProjAngles)
vertical = (fabsf(angles[a].fRayX) <= fabsf(angles[a].fRayY));
if (a == dims.iProjAngles || vertical != blockVertical) {
// block done
blockEnd = a;
if (blockStart != blockEnd) {
dim3 dimGrid((blockEnd-blockStart+g_anglesPerBlock-1)/g_anglesPerBlock,
(dims.iProjDets+g_detBlockSize-1)/g_detBlockSize); // angle blocks, detector blocks
// TODO: check if we can't immediately
// destroy the stream after use
cudaStream_t stream;
cudaStreamCreate(&stream);
streams.push_back(stream);
//printf("angle block: %d to %d, %d\n", blockStart, blockEnd, blockVertical);
if (!blockVertical)
for (unsigned int i = 0; i < dims.iVolWidth; i += g_blockSlices)
FPhorizontal_simple<<<dimGrid, dimBlock, 0, stream>>>(D_projData, projPitch, i, blockStart, blockEnd, dims, outputScale);
else
for (unsigned int i = 0; i < dims.iVolHeight; i += g_blockSlices)
FPvertical_simple<<<dimGrid, dimBlock, 0, stream>>>(D_projData, projPitch, i, blockStart, blockEnd, dims, outputScale);
}
blockVertical = vertical;
blockStart = a;
}
}
for (std::list<cudaStream_t>::iterator iter = streams.begin(); iter != streams.end(); ++iter)
cudaStreamDestroy(*iter);
streams.clear();
cudaThreadSynchronize();
cudaTextForceKernelsCompletion();
cudaFreeArray(D_dataArray);
return true;
}
bool FP_simple(float* D_volumeData, unsigned int volumePitch,
float* D_projData, unsigned int projPitch,
const SDimensions& dims, const SParProjection* angles,
float outputScale)
{
for (unsigned int iAngle = 0; iAngle < dims.iProjAngles; iAngle += g_MaxAngles) {
SDimensions subdims = dims;
unsigned int iEndAngle = iAngle + g_MaxAngles;
if (iEndAngle >= dims.iProjAngles)
iEndAngle = dims.iProjAngles;
subdims.iProjAngles = iEndAngle - iAngle;
bool ret;
ret = FP_simple_internal(D_volumeData, volumePitch,
D_projData + iAngle * projPitch, projPitch,
subdims, angles + iAngle,
outputScale);
if (!ret)
return false;
}
return true;
}
bool FP(float* D_volumeData, unsigned int volumePitch,
float* D_projData, unsigned int projPitch,
const SDimensions& dims, const SParProjection* angles,
float outputScale)
{
return FP_simple(D_volumeData, volumePitch, D_projData, projPitch,
dims, angles, outputScale);
}
}
#ifdef STANDALONE
using namespace astraCUDA;
int main()
{
float* D_volumeData;
float* D_projData;
SDimensions dims;
dims.iVolWidth = 1024;
dims.iVolHeight = 1024;
dims.iProjAngles = 512;
dims.iProjDets = 1536;
dims.fDetScale = 1.0f;
dims.iRaysPerDet = 1;
unsigned int volumePitch, projPitch;
allocateVolume(D_volumeData, dims.iVolWidth, dims.iVolHeight, volumePitch);
printf("pitch: %u\n", volumePitch);
allocateVolume(D_projData, dims.iProjDets, dims.iProjAngles, projPitch);
printf("pitch: %u\n", projPitch);
unsigned int y, x;
float* img = loadImage("phantom.png", y, x);
float* sino = new float[dims.iProjAngles * dims.iProjDets];
memset(sino, 0, dims.iProjAngles * dims.iProjDets * sizeof(float));
copyVolumeToDevice(img, dims.iVolWidth, dims.iVolWidth, dims.iVolHeight, D_volumeData, volumePitch);
copySinogramToDevice(sino, dims.iProjDets, dims.iProjDets, dims.iProjAngles, D_projData, projPitch);
float* angle = new float[dims.iProjAngles];
for (unsigned int i = 0; i < dims.iProjAngles; ++i)
angle[i] = i*(M_PI/dims.iProjAngles);
FP(D_volumeData, volumePitch, D_projData, projPitch, dims, angle, 0, 1.0f);
delete[] angle;
copySinogramFromDevice(sino, dims.iProjDets, dims.iProjDets, dims.iProjAngles, D_projData, projPitch);
float s = 0.0f;
for (unsigned int y = 0; y < dims.iProjAngles; ++y)
for (unsigned int x = 0; x < dims.iProjDets; ++x)
s += sino[y*dims.iProjDets+x] * sino[y*dims.iProjDets+x];
printf("cpu norm: %f\n", s);
//zeroVolume(D_projData, projPitch, dims.iProjDets, dims.iProjAngles);
s = dotProduct2D(D_projData, projPitch, dims.iProjDets, dims.iProjAngles, 1, 0);
printf("gpu norm: %f\n", s);
saveImage("sino.png",dims.iProjAngles,dims.iProjDets,sino);
return 0;
}
#endif
|