summaryrefslogtreecommitdiffstats
path: root/cuda/2d/sart.cu
blob: aff77a31bd4df3982c178c9b328c0e574fdd999d (plain)
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
/*
-----------------------------------------------------------------------
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/sart.h"
#include "astra/cuda/2d/util.h"
#include "astra/cuda/2d/arith.h"
#include "astra/cuda/2d/fan_fp.h"
#include "astra/cuda/2d/fan_bp.h"
#include "astra/cuda/2d/par_fp.h"
#include "astra/cuda/2d/par_bp.h"

#include <cstdio>
#include <cassert>

namespace astraCUDA {

// FIXME: Remove these functions. (Outdated)
__global__ void devMUL_SART(float* pfOut, const float* pfIn, unsigned int pitch, unsigned int width)
{
	unsigned int x = threadIdx.x + 16*blockIdx.x;
	if (x >= width) return;

	pfOut[x] *= pfIn[x];
}

void MUL_SART(float* pfOut, const float* pfIn, unsigned int pitch, unsigned int width)
{
	dim3 blockSize(16,16);
	dim3 gridSize((width+15)/16, 1);

	devMUL_SART<<<gridSize, blockSize>>>(pfOut, pfIn, pitch, width);

	cudaTextForceKernelsCompletion();
}



SART::SART() : ReconAlgo()
{
	D_projData = 0;
	D_tmpData = 0;

	D_lineWeight = 0;

	projectionOrder = 0;
	projectionCount = 0;
	iteration = 0;
	customOrder = false;

	fRelaxation = 1.0f;
}


SART::~SART()
{
	reset();
}

void SART::reset()
{
	cudaFree(D_projData);
	cudaFree(D_tmpData);
	cudaFree(D_lineWeight);

	D_projData = 0;
	D_tmpData = 0;

	D_lineWeight = 0;

	useVolumeMask = false;
	useSinogramMask = false;

	if (projectionOrder != NULL) delete[] projectionOrder;
	projectionOrder = 0;
	projectionCount = 0;
	iteration = 0;
	customOrder = false;
	fRelaxation = 1.0f;

	ReconAlgo::reset();
}

bool SART::init()
{
	if (useVolumeMask) {
		allocateVolumeData(D_tmpData, tmpPitch, dims);
		zeroVolumeData(D_tmpData, tmpPitch, dims);
	}

	// NB: Non-standard dimensions
	SDimensions linedims = dims;
	linedims.iProjAngles = 1;
	allocateProjectionData(D_projData, projPitch, linedims);
	zeroProjectionData(D_projData, projPitch, linedims);
	
	allocateProjectionData(D_lineWeight, linePitch, dims);
	zeroProjectionData(D_lineWeight, linePitch, dims);

	// We can't precompute lineWeights when using a mask
	if (!useVolumeMask)
		precomputeWeights();

	// TODO: check if allocations succeeded
	return true;
}

bool SART::setProjectionOrder(int* _projectionOrder, int _projectionCount)
{
	customOrder = true;
	projectionCount = _projectionCount;
	projectionOrder = new int[projectionCount];
	for (int i = 0; i < projectionCount; i++) {
		projectionOrder[i] = _projectionOrder[i];
	}

	return true;
}


bool SART::precomputeWeights()
{
	zeroProjectionData(D_lineWeight, linePitch, dims);
	if (useVolumeMask) {
		callFP(D_maskData, maskPitch, D_lineWeight, linePitch, 1.0f);
	} else {
		// Allocate tmpData temporarily
		allocateVolumeData(D_tmpData, tmpPitch, dims);
		zeroVolumeData(D_tmpData, tmpPitch, dims);


		processVol<opSet>(D_tmpData, 1.0f, tmpPitch, dims);
		callFP(D_tmpData, tmpPitch, D_lineWeight, linePitch, 1.0f);


		cudaFree(D_tmpData);
		D_tmpData = 0;
	}
	processSino<opInvert>(D_lineWeight, linePitch, dims);

	return true;
}

bool SART::iterate(unsigned int iterations)
{
	if (useVolumeMask)
		precomputeWeights();

	// iteration
	for (unsigned int iter = 0; iter < iterations && !astra::shouldAbort(); ++iter) {

		int angle;
		if (customOrder) {
			angle = projectionOrder[iteration % projectionCount];
		} else {
			angle = iteration % dims.iProjAngles;  
		}

		// copy one line of sinogram to projection data
		// NB: Non-standard dimensions
		SDimensions linedims = dims;
		linedims.iProjAngles = 1;
		duplicateProjectionData(D_projData, D_sinoData + angle*sinoPitch, sinoPitch, linedims);

		// do FP, subtracting projection from sinogram
		if (useVolumeMask) {
				duplicateVolumeData(D_tmpData, D_volumeData, volumePitch, dims);
				processVol<opMul>(D_tmpData, D_maskData, tmpPitch, dims);
				callFP_SART(D_tmpData, tmpPitch, D_projData, projPitch, angle, -1.0f);
		} else {
				callFP_SART(D_volumeData, volumePitch, D_projData, projPitch, angle, -1.0f);
		}

		MUL_SART(D_projData, D_lineWeight + angle*linePitch, projPitch, dims.iProjDets);

		if (useVolumeMask) {
			// BP, mask, and add back
			// TODO: Try putting the masking directly in the BP
			zeroVolumeData(D_tmpData, tmpPitch, dims);
			callBP_SART(D_tmpData, tmpPitch, D_projData, projPitch, angle, fRelaxation);
			processVol<opAddMul>(D_volumeData, D_maskData, D_tmpData, volumePitch, dims);
		} else {
			callBP_SART(D_volumeData, volumePitch, D_projData, projPitch, angle, fRelaxation);
		}

		if (useMinConstraint)
			processVol<opClampMin>(D_volumeData, fMinConstraint, volumePitch, dims);
		if (useMaxConstraint)
			processVol<opClampMax>(D_volumeData, fMaxConstraint, volumePitch, dims);

		iteration++;

	}

	return true;
}

float SART::computeDiffNorm()
{
	unsigned int pPitch;
	float *D_p;
	allocateProjectionData(D_p, pPitch, dims);

	// copy sinogram to D_p
	duplicateProjectionData(D_p, D_sinoData, sinoPitch, dims);

	// do FP, subtracting projection from sinogram
	if (useVolumeMask) {
			duplicateVolumeData(D_tmpData, D_volumeData, volumePitch, dims);
			processVol<opMul>(D_tmpData, D_maskData, tmpPitch, dims);
			callFP(D_tmpData, tmpPitch, D_p, pPitch, -1.0f);
	} else {
			callFP(D_volumeData, volumePitch, D_p, pPitch, -1.0f);
	}


	// compute norm of D_p
	float s = dotProduct2D(D_p, pPitch, dims.iProjDets, dims.iProjAngles);

	cudaFree(D_p);

	return sqrt(s);
}

bool SART::callFP_SART(float* D_volumeData, unsigned int volumePitch,
                       float* D_projData, unsigned int projPitch,
                       unsigned int angle, float outputScale)
{
	SDimensions d = dims;
	d.iProjAngles = 1;
	if (parProjs) {
		assert(!fanProjs);
		return FP(D_volumeData, volumePitch, D_projData, projPitch,
		          d, &parProjs[angle], outputScale * fProjectorScale);
	} else {
		assert(fanProjs);
		return FanFP(D_volumeData, volumePitch, D_projData, projPitch,
		             d, &fanProjs[angle], outputScale * fProjectorScale);
	}
}

bool SART::callBP_SART(float* D_volumeData, unsigned int volumePitch,
                       float* D_projData, unsigned int projPitch,
                       unsigned int angle, float outputScale)
{
	if (parProjs) {
		assert(!fanProjs);
		return BP_SART(D_volumeData, volumePitch, D_projData, projPitch,
		               angle, dims, parProjs, outputScale * fProjectorScale);
	} else {
		assert(fanProjs);
		return FanBP_SART(D_volumeData, volumePitch, D_projData, projPitch,
		                  angle, dims, fanProjs, outputScale * fProjectorScale);
	}

}


}