bioem_cuda.cu 14.1 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
#define BIOEM_GPUCODE

#if defined(_WIN32)
#include <windows.h>
#endif

#include <iostream>
using namespace std;

#include "bioem_cuda_internal.h"
#include "bioem_algorithm.h"

13
14
15
16
17
18
19
#define checkCudaErrors(error) \
{ \
	if ((error) != cudaSuccess) \
	{ \
		printf("CUDA Error %d / %s (%s: %d)\n", error, cudaGetErrorString(error), __FILE__, __LINE__); \
		exit(1); \
	} \
20
21
}

David Rohr's avatar
David Rohr committed
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
static const char *cufftGetErrorStrung(cufftResult error)
{
    switch (error)
    {
        case CUFFT_SUCCESS:
            return "CUFFT_SUCCESS";

        case CUFFT_INVALID_PLAN:
            return "CUFFT_INVALID_PLAN";

        case CUFFT_ALLOC_FAILED:
            return "CUFFT_ALLOC_FAILED";

        case CUFFT_INVALID_TYPE:
            return "CUFFT_INVALID_TYPE";

        case CUFFT_INVALID_VALUE:
            return "CUFFT_INVALID_VALUE";

        case CUFFT_INTERNAL_ERROR:
            return "CUFFT_INTERNAL_ERROR";

        case CUFFT_EXEC_FAILED:
            return "CUFFT_EXEC_FAILED";

        case CUFFT_SETUP_FAILED:
            return "CUFFT_SETUP_FAILED";

        case CUFFT_INVALID_SIZE:
            return "CUFFT_INVALID_SIZE";

        case CUFFT_UNALIGNED_DATA:
            return "CUFFT_UNALIGNED_DATA";
    }
    return "UNKNOWN";
}

59
60
61
62
bioem_cuda::bioem_cuda()
{
	deviceInitialized = 0;
	GPUAlgo = getenv("GPUALGO") == NULL ? 2 : atoi(getenv("GPUALGO"));
63
64
	GPUAsync = getenv("GPUASYNC") == NULL ? 1 : atoi(getenv("GPUASYNC"));
	GPUWorkload = getenv("GPUWORKLOAD") == NULL ? 100 : atoi(getenv("GPUWORKLOAD"));
65
	GPUDualStream = getenv("GPUDUALSTREAM") == NULL ? 1 : atoi(getenv("GPUDUALSTREAM"));
66
67
68
69
70
71
72
}

bioem_cuda::~bioem_cuda()
{
	deviceExit();
}

73
__global__ void compareRefMap_kernel(const int iOrient, const int iConv, const myfloat_t* pMap, bioem_Probability pProb, const bioem_param_device param, const bioem_RefMap_Mod RefMap, const int cent_x, const int cent_y, const int maxRef)
74
75
{
	const int iRefMap = myBlockIdxX * myBlockDimX + myThreadIdxX;
76
	if (iRefMap < maxRef)
77
	{
78
		compareRefMap<0>(iRefMap, iOrient, iConv, pMap, pProb, param, RefMap, cent_x, cent_y);
79
80
81
	}
}

82
__global__ void compareRefMapShifted_kernel(const int iOrient, const int iConv, const myfloat_t* pMap, bioem_Probability pProb, const bioem_param_device param, const bioem_RefMap_Mod RefMap, const int maxRef)
83
84
{
	const int iRefMap = myBlockIdxX * myBlockDimX + myThreadIdxX;
85
	if (iRefMap < maxRef)
86
	{
87
		compareRefMapShifted<1>(iRefMap, iOrient, iConv, pMap, pProb, param, RefMap);
88
89
90
	}
}

91
92
93
94
95
96
__global__ void cudaZeroMem(void* ptr, size_t size)
{
	int* myptr = (int*) ptr;
	int mysize = size / sizeof(int);
	int myid = myBlockDimX * myBlockIdxX + myThreadIdxX;
	int mygrid = myBlockDimX * myGridDimX;
97
	for (int i = myid; i < mysize; i += mygrid) myptr[i] = 0;
98
99
}

100
__global__ void compareRefMapLoopShifts_kernel(const int iOrient, const int iConv, const myfloat_t* pMap, bioem_Probability pProb, const bioem_param_device param, const bioem_RefMap RefMap, const int blockoffset, const int nShifts, const int nShiftBits, const int maxRef)
101
102
103
104
105
106
107
108
109
{
	const size_t myid = (myBlockIdxX + blockoffset) * myBlockDimX + myThreadIdxX;
	const int iRefMap = myid >> (nShiftBits << 1);
	const int myRef = myThreadIdxX >> (nShiftBits << 1);
	const int myShiftIdx = (myid >> nShiftBits) & (nShifts - 1);
	const int myShiftIdy = myid & (nShifts - 1);
	const int myShift = myid & (nShifts * nShifts - 1);
	const int cent_x = myShiftIdx * param.GridSpaceCenter - param.maxDisplaceCenter;
	const int cent_y = myShiftIdy * param.GridSpaceCenter - param.maxDisplaceCenter;
110

111
	const bool threadActive = myShiftIdx < nShifts && myShiftIdy < nShifts && iRefMap < maxRef;
112

113
	compareRefMap<2>(iRefMap, iOrient, iConv, pMap, pProb, param, RefMap, cent_x, cent_y, myShift, nShifts * nShifts, myRef, threadActive);
114
115
}

116
__global__ void multComplexMap(const mycomplex_t* convmap, const mycomplex_t* refmap, mycuComplex_t* out, const int NumberPixelsTotal, const int MapSize, const int NumberMaps, const int Offset)
117
118
{
	if (myBlockIdxX >= NumberMaps) return;
119
	const mycuComplex_t* myin = (mycuComplex_t*) &refmap[(myBlockIdxX + Offset) * MapSize];
120
	const mycuComplex_t* myconv = (mycuComplex_t*) convmap;
121
	mycuComplex_t* myout = &out[myBlockIdxX * MapSize];
122
	for(int i = myThreadIdxX; i < NumberPixelsTotal; i += myBlockDimX)
123
	{
124
125
126
127
128
129
130
		mycuComplex_t val;
		const mycuComplex_t conv = myconv[i];
		const mycuComplex_t in = myin[i];

		val.x = conv.x * in.x + conv.y * in.y;
		val.y = conv.y * in.x - conv.x * in.y;
		myout[i] = val;
131
132
133
	}
}

134
__global__ void cuDoRefMapsFFT(const int iOrient, const int iConv, const myfloat_t* lCC, const myfloat_t sumC, const myfloat_t sumsquareC, bioem_Probability pProb, const bioem_param_device param, const bioem_RefMap RefMap, const int maxRef, const int Offset)
135
{
136
	if (myBlockIdxX * myBlockDimX + myThreadIdxX >= maxRef) return;
137
138
	const int iRefMap = myBlockIdxX * myBlockDimX + myThreadIdxX + Offset;
	const myfloat_t* mylCC = &lCC[(myBlockIdxX * myBlockDimX + myThreadIdxX) * param.NumberPixels * param.NumberPixels];
139
	doRefMapFFT(iRefMap, iOrient, iConv, mylCC, sumC, sumsquareC, pProb, param, RefMap);
140
141
}

142
143
144
145
146
147
148
149
150
151
152
153
154
template <class T> static inline T divup(T num, T divider) {return((num + divider - 1) / divider);}
static inline bool IsPowerOf2(int x) {return ((x > 0) && ((x & (x - 1)) == 0));}
#if defined(_WIN32)
static inline int ilog2 (int value)
{
	DWORD index;
	_BitScanReverse (&index, value);
	return(value);
}
#else
static inline int ilog2(int value) {return 31 - __builtin_clz(value);}
#endif

155
int bioem_cuda::compareRefMaps(int iOrient, int iConv, const myfloat_t* conv_map, mycomplex_t* localmultFFT, myfloat_t sumC, myfloat_t sumsquareC, const int startMap)
156
{
157
158
159
160
161
	if (startMap)
	{
		cout << "Error startMap not implemented for GPU Code\n";
		exit(1);
	}
162
163
164
165
	if (GPUAsync)
	{
		checkCudaErrors(cudaEventSynchronize(cudaEvent[iConv & 1]));
	}
166

167
	if (FFTAlgo)
168
	{
169
170
		checkCudaErrors(cudaMemcpyAsync(&pConvMapFFT[(iConv & 1) * param.FFTMapSize], localmultFFT, param.FFTMapSize * sizeof(mycomplex_t), cudaMemcpyHostToDevice, cudaStream[0]));
		if (GPUDualStream)
171
		{
172
173
174
175
176
177
			checkCudaErrors(cudaEventRecord(cudaFFTEvent[0], cudaStream[0]));
			checkCudaErrors(cudaStreamWaitEvent(cudaStream[1], cudaFFTEvent[0], 0));
		}
		for (int i = 0, j = 0; i < maxRef; i += CUDA_FFTS_AT_ONCE, j++)
		{
			if (!GPUDualStream) j = 0;
178
			const int num = min(CUDA_FFTS_AT_ONCE, maxRef - i);
179
180
			multComplexMap<<<num, CUDA_THREAD_COUNT, 0, cudaStream[j & 1]>>>(&pConvMapFFT[(iConv & 1) * param.FFTMapSize], pRefMapsFFT, pFFTtmp2[j & 1], param.param_device.NumberPixels * param.param_device.NumberFFTPixels1D, param.FFTMapSize, num, i);
			cufftResult err = mycufftExecC2R(i + CUDA_FFTS_AT_ONCE > maxRef ? plan[1][j & 1] : plan[0][j & 1], pFFTtmp2[j & 1], pFFTtmp[j & 1]);
David Rohr's avatar
David Rohr committed
181
			if (err != CUFFT_SUCCESS)
182
			{
David Rohr's avatar
David Rohr committed
183
				cout << "Error running CUFFT " << cufftGetErrorStrung(err) << "\n";
184
185
				exit(1);
			}
186
			cuDoRefMapsFFT<<<divup(num, CUDA_THREAD_COUNT), CUDA_THREAD_COUNT, 0, cudaStream[j & 1]>>>(iOrient, iConv, pFFTtmp[j & 1], sumC, sumsquareC, pProb_device, param.param_device, *gpumap, num, i);
187
		}
188
		checkCudaErrors(cudaGetLastError());
189
190
191
192
193
		if (GPUDualStream)
		{
			checkCudaErrors(cudaEventRecord(cudaFFTEvent[1], cudaStream[1]));
			checkCudaErrors(cudaStreamWaitEvent(cudaStream[0], cudaFFTEvent[1], 0));
		}
194
195
196
	}
	else
	{
197
		checkCudaErrors(cudaMemcpyAsync(pConvMap_device[iConv & 1], conv_map, sizeof(myfloat_t) * RefMap.refMapSize, cudaMemcpyHostToDevice, cudaStream[0]));
198
199

		if (GPUAlgo == 2) //Loop over shifts
200
		{
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
			const int nShifts = 2 * param.param_device.maxDisplaceCenter / param.param_device.GridSpaceCenter + 1;
			if (!IsPowerOf2(nShifts))
			{
				cout << "Invalid number of displacements, no power of two\n";
				exit(1);
			}
			if (CUDA_THREAD_COUNT % (nShifts * nShifts))
			{
				cout << "CUDA Thread count (" << CUDA_THREAD_COUNT << ") is no multiple of number of shifts (" << (nShifts * nShifts) << ")\n";
				exit(1);
			}
			if (nShifts > CUDA_MAX_SHIFT_REDUCE)
			{
				cout << "Too many displacements for CUDA reduction\n";
				exit(1);
			}
			const int nShiftBits = ilog2(nShifts);
			size_t totalBlocks = divup((size_t) maxRef * (size_t) nShifts * (size_t) nShifts, (size_t) CUDA_THREAD_COUNT);
			size_t nBlocks = CUDA_BLOCK_COUNT;
220
			for (size_t i = 0; i < totalBlocks; i += nBlocks)
221
			{
222
				compareRefMapLoopShifts_kernel<<<min(nBlocks, totalBlocks - i), CUDA_THREAD_COUNT, (CUDA_THREAD_COUNT * 2 + CUDA_THREAD_COUNT / (nShifts * nShifts) * 4) * sizeof(myfloat_t), cudaStream[0] >>> (iOrient, iConv, pConvMap_device[iConv & 1], pProb_device, param.param_device, *gpumap, i, nShifts, nShiftBits, maxRef);
223
			}
224
		}
225
		else if (GPUAlgo == 1) //Split shifts in multiple kernels
226
		{
227
			for (int cent_x = -param.param_device.maxDisplaceCenter; cent_x <= param.param_device.maxDisplaceCenter; cent_x = cent_x + param.param_device.GridSpaceCenter)
228
			{
229
				for (int cent_y = -param.param_device.maxDisplaceCenter; cent_y <= param.param_device.maxDisplaceCenter; cent_y = cent_y + param.param_device.GridSpaceCenter)
230
				{
231
					compareRefMap_kernel<<<divup(maxRef, CUDA_THREAD_COUNT), CUDA_THREAD_COUNT, 0, cudaStream[0]>>> (iOrient, iConv, pConvMap_device[iConv & 1], pProb_device, param.param_device, *pRefMap_device_Mod, cent_x, cent_y, maxRef);
232
233
				}
			}
234
		}
235
		else if (GPUAlgo == 0) //All shifts in one kernel
236
		{
237
			compareRefMapShifted_kernel<<<divup(maxRef, CUDA_THREAD_COUNT), CUDA_THREAD_COUNT, 0, cudaStream[0]>>> (iOrient, iConv, pConvMap_device[iConv & 1], pProb_device, param.param_device, *pRefMap_device_Mod, maxRef);
238
		}
239
		else
240
		{
241
242
			cout << "Invalid GPU Algorithm selected\n";
			exit(1);
243
		}
244
	}
245
246
	if (GPUWorkload < 100)
	{
247
		bioem::compareRefMaps(iOrient, iConv, conv_map, localmultFFT, sumC, sumsquareC, maxRef);
248
	}
249
250
	if (GPUAsync)
	{
251
		checkCudaErrors(cudaEventRecord(cudaEvent[iConv & 1], cudaStream[0]));
252
	}
253
254
	else
	{
255
		checkCudaErrors(cudaStreamSynchronize(cudaStream[0]));
256
257
258
259
260
261
262
	}
	return(0);
}

int bioem_cuda::deviceInit()
{
	deviceExit();
263

264
265
	if (FFTAlgo) GPUAlgo = 2;

266
267
268
269
270
	gpumap = new bioem_RefMap;
	memcpy(gpumap, &RefMap, sizeof(bioem_RefMap));
	if (FFTAlgo == 0)
	{
		checkCudaErrors(cudaMalloc(&maps, sizeof(myfloat_t) * RefMap.ntotRefMap * RefMap.refMapSize));
271
272
273
274
275
276
277
278
279
280
281
282
283

		if (GPUAlgo == 0 || GPUAlgo == 1)
		{
			pRefMap_device_Mod = (bioem_RefMap_Mod*) gpumap;
			bioem_RefMap_Mod* RefMapGPU = new bioem_RefMap_Mod;
			RefMapGPU->init(RefMap);
			checkCudaErrors(cudaMemcpy(maps, RefMapGPU->maps, sizeof(myfloat_t) * RefMap.ntotRefMap * RefMap.refMapSize, cudaMemcpyHostToDevice));
			delete RefMapGPU;
		}
		else
		{
			checkCudaErrors(cudaMemcpy(maps, RefMap.maps, sizeof(myfloat_t) * RefMap.ntotRefMap * RefMap.refMapSize, cudaMemcpyHostToDevice));
		}
284
285
286
287
288
289
290
291
292
	}
	checkCudaErrors(cudaMalloc(&sum, sizeof(myfloat_t) * RefMap.ntotRefMap));
	checkCudaErrors(cudaMemcpy(sum, RefMap.sum_RefMap, sizeof(myfloat_t) * RefMap.ntotRefMap, cudaMemcpyHostToDevice));
	checkCudaErrors(cudaMalloc(&sumsquare, sizeof(myfloat_t) * RefMap.ntotRefMap));
	checkCudaErrors(cudaMemcpy(sumsquare, RefMap.sumsquare_RefMap, sizeof(myfloat_t) * RefMap.ntotRefMap, cudaMemcpyHostToDevice));
	gpumap->maps = maps;
	gpumap->sum_RefMap = sum;
	gpumap->sumsquare_RefMap = sumsquare;

293
294
295
	pProb_device = pProb;
	checkCudaErrors(cudaMalloc(&pProb_device.ptr, pProb_device.get_size(RefMap.ntotRefMap, param.nTotGridAngles)));
	pProb_device.set_pointers();
296
	for (int i = 0; i < 2; i++)
297
	{
298
		checkCudaErrors(cudaStreamCreate(&cudaStream[i]));
299
		checkCudaErrors(cudaEventCreate(&cudaEvent[i]));
300
		checkCudaErrors(cudaEventCreate(&cudaFFTEvent[i]));
301
		checkCudaErrors(cudaMalloc(&pConvMap_device[i], sizeof(myfloat_t) * RefMap.refMapSize));
302
	}
303

304
305
	if (FFTAlgo)
	{
306
		checkCudaErrors(cudaMalloc(&pRefMapsFFT, RefMap.ntotRefMap * param.FFTMapSize * sizeof(mycomplex_t)));
307
308
309
310
		checkCudaErrors(cudaMalloc(&pFFTtmp2[0], CUDA_FFTS_AT_ONCE * param.FFTMapSize * 2 * sizeof(mycomplex_t)));
		checkCudaErrors(cudaMalloc(&pFFTtmp[0], CUDA_FFTS_AT_ONCE * param.param_device.NumberPixels * param.param_device.NumberPixels * 2 * sizeof(myfloat_t)));
		pFFTtmp2[1] = pFFTtmp2[0] + CUDA_FFTS_AT_ONCE * param.FFTMapSize;
		pFFTtmp[1] = pFFTtmp[0] + CUDA_FFTS_AT_ONCE * param.param_device.NumberPixels * param.param_device.NumberPixels;
311
312
		checkCudaErrors(cudaMalloc(&pConvMapFFT, param.FFTMapSize * sizeof(mycomplex_t) * 2));
		checkCudaErrors(cudaMemcpy(pRefMapsFFT, RefMap.RefMapsFFT, RefMap.ntotRefMap * param.FFTMapSize * sizeof(mycomplex_t), cudaMemcpyHostToDevice));
313
314
	}

315
316
317
318
319
320
321
	deviceInitialized = 1;
	return(0);
}

int bioem_cuda::deviceExit()
{
	if (deviceInitialized == 0) return(0);
322

323

324
	cudaFree(pProb_device.ptr);
325
326
	cudaFree(sum);
	cudaFree(sumsquare);
327
	for (int i = 0; i < 2; i++)
328
	{
329
		cudaStreamDestroy(cudaStream[i]);
330
		cudaEventDestroy(cudaEvent[i]);
331
		cudaEventDestroy(cudaFFTEvent[i]);
332
		cudaFree(pConvMap_device[i]);
333
	}
334
335
336
337
	if (FFTAlgo)
	{
		cudaFree(pRefMapsFFT);
		cudaFree(pConvMapFFT);
338
339
		cudaFree(pFFTtmp[0]);
		cudaFree(pFFTtmp2[0]);
340
	}
341
342
343
344
345
346
347
348
349
	else
	{
		cudaFree(maps);
	}
	if (GPUAlgo == 0 || GPUAlgo == 1)
	{
		cudaFree(pRefMap_device_Mod);
	}
	delete gpumap;
350
	cudaThreadExit();
351

352
353
354
355
356
357
	deviceInitialized = 0;
	return(0);
}

int bioem_cuda::deviceStartRun()
{
358
	maxRef = GPUWorkload >= 100 ? RefMap.ntotRefMap : ((size_t) RefMap.ntotRefMap * (size_t) GPUWorkload / 100);
359

360
	cudaMemcpy(pProb_device.ptr, pProb.ptr, pProb.get_size(RefMap.ntotRefMap, param.nTotGridAngles), cudaMemcpyHostToDevice);
361
362
363

	if (FFTAlgo)
	{
364
		for (int j = 0;j < 2;j++)
365
		{
366
			for (int i = 0; i < 2; i++)
367
			{
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
				int n[2] = {param.param_device.NumberPixels, param.param_device.NumberPixels};
				if (cufftPlanMany(&plan[i][j], 2, n, NULL, 1, param.FFTMapSize, NULL, 1, 0, MY_CUFFT_C2R, i ? (maxRef % CUDA_FFTS_AT_ONCE) : CUDA_FFTS_AT_ONCE) != CUFFT_SUCCESS)
				{
					cout << "Error planning CUFFT\n";
					exit(1);
				}
				if (cufftSetCompatibilityMode(plan[i][j], CUFFT_COMPATIBILITY_NATIVE) != CUFFT_SUCCESS)
				{
					cout << "Error planning CUFFT compatibility\n";
					exit(1);
				}
				if (cufftSetStream(plan[i][j], cudaStream[j]) != CUFFT_SUCCESS)
				{
					cout << "Error setting CUFFT stream\n";
					exit(1);
				}
384
			}
385
			if (!GPUDualStream) break;
386
387
		}
	}
388
389
390
391
392
	return(0);
}

int bioem_cuda::deviceFinishRun()
{
393
	if (GPUAsync) cudaStreamSynchronize(cudaStream[0]);
394
	cudaMemcpy(pProb.ptr, pProb_device.ptr, pProb.get_size(RefMap.ntotRefMap, param.nTotGridAngles), cudaMemcpyDeviceToHost);
395

396
397
	if (FFTAlgo)
	{
398
399
400
401
402
		for (int j = 0;j < 2;j++)
		{
			for (int i = 0; i < 2; i++) cufftDestroy(plan[i][j]);
			if (!GPUDualStream) break;
		}
403
404
	}

405
406
407
408
409
410
411
	return(0);
}

bioem* bioem_cuda_create()
{
	return new bioem_cuda;
}