bioem_cuda.cu 13 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 66 67 68 69 70 71
}

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

72
__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)
73 74
{
	const int iRefMap = myBlockIdxX * myBlockDimX + myThreadIdxX;
75
	if (iRefMap < maxRef)
76
	{
77
		compareRefMap<0>(iRefMap, iOrient, iConv, pMap, pProb, param, RefMap, cent_x, cent_y);
78 79 80
	}
}

81
__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)
82 83
{
	const int iRefMap = myBlockIdxX * myBlockDimX + myThreadIdxX;
84
	if (iRefMap < maxRef)
85
	{
86
		compareRefMapShifted<1>(iRefMap, iOrient, iConv, pMap, pProb, param, RefMap);
87 88 89
	}
}

90 91 92 93 94 95
__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;
96
	for (int i = myid; i < mysize; i += mygrid) myptr[i] = 0;
97 98
}

99
__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)
100 101 102 103 104 105 106 107 108
{
	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;
109

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

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

115
__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)
116 117
{
	if (myBlockIdxX >= NumberMaps) return;
118
	const mycomplex_t* myin = &refmap[myBlockIdxX * MapSize + Offset];
119
	mycuComplex_t* myout = &out[myBlockIdxX * MapSize];
120
	for(int i = myThreadIdxX; i < NumberPixelsTotal; i += myBlockDimX)
121 122 123 124 125 126
	{
		myout[i].x = convmap[i][0] * myin[i][0] + convmap[i][1] * myin[i][1];
		myout[i].y = convmap[i][1] * myin[i][0] - convmap[i][0] * myin[i][1];
	}
}

127
__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)
128
{
129 130
	const int iRefMap = myBlockIdxX * myBlockDimX + myThreadIdxX + Offset;
	const myfloat_t* mylCC = &lCC[(myBlockIdxX * myBlockDimX + myThreadIdxX) * param.NumberPixels * param.NumberPixels];
131
	if (iRefMap >= maxRef) return;
132
	doRefMapFFT(iRefMap, iOrient, iConv, mylCC, sumC, sumsquareC, pProb, param, RefMap);
133 134
}

135 136 137 138 139 140 141 142 143 144 145 146 147
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

148
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)
149
{
150 151 152 153 154
	if (startMap)
	{
		cout << "Error startMap not implemented for GPU Code\n";
		exit(1);
	}
155 156 157 158
	if (GPUAsync)
	{
		checkCudaErrors(cudaEventSynchronize(cudaEvent[iConv & 1]));
	}
159

160
	if (FFTAlgo)
161
	{
162
		checkCudaErrors(cudaMemcpyAsync(&pConvMapFFT[(iConv & 1) * param.FFTMapSize], localmultFFT, param.FFTMapSize * sizeof(mycomplex_t), cudaMemcpyHostToDevice, cudaStream));
163
		for (int i = 0; i < maxRef; i += CUDA_FFTS_AT_ONCE)
164
		{
165
			const int num = min(CUDA_FFTS_AT_ONCE, maxRef - i);
166
			multComplexMap<<<num, CUDA_THREAD_COUNT, 0, cudaStream>>>(&pConvMapFFT[(iConv & 1) * param.FFTMapSize], pRefMapsFFT, pFFTtmp2, param.param_device.NumberPixels * param.param_device.NumberFFTPixels1D, param.FFTMapSize, num, i);
David Rohr's avatar
David Rohr committed
167 168
			cufftResult err = mycufftExecC2R(i + CUDA_FFTS_AT_ONCE > maxRef ? plan[1] : plan[0], pFFTtmp2, pFFTtmp);
			if (err != CUFFT_SUCCESS)
169
			{
David Rohr's avatar
David Rohr committed
170
				cout << "Error running CUFFT " << cufftGetErrorStrung(err) << "\n";
171 172
				exit(1);
			}
173
			cuDoRefMapsFFT<<<divup(num, CUDA_THREAD_COUNT), CUDA_THREAD_COUNT, 0, cudaStream>>>(iOrient, iConv, pFFTtmp, sumC, sumsquareC, pProb_device, param.param_device, *gpumap, num, i);
174
		}
175 176 177 178
		checkCudaErrors(cudaGetLastError());
	}
	else
	{
179
		checkCudaErrors(cudaMemcpyAsync(pConvMap_device[iConv & 1], conv_map, sizeof(myfloat_t) * RefMap.refMapSize, cudaMemcpyHostToDevice, cudaStream));
180 181

		if (GPUAlgo == 2) //Loop over shifts
182
		{
183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201
			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;
202
			for (size_t i = 0; i < totalBlocks; i += nBlocks)
203
			{
204
				compareRefMapLoopShifts_kernel<<<min(nBlocks, totalBlocks - i), CUDA_THREAD_COUNT, (CUDA_THREAD_COUNT * 2 + CUDA_THREAD_COUNT / (nShifts * nShifts) * 4) * sizeof(myfloat_t), cudaStream >>> (iOrient, iConv, pConvMap_device[iConv & 1], pProb_device, param.param_device, *gpumap, i, nShifts, nShiftBits, maxRef);
205
			}
206
		}
207
		else if (GPUAlgo == 1) //Split shifts in multiple kernels
208
		{
209
			for (int cent_x = -param.param_device.maxDisplaceCenter; cent_x <= param.param_device.maxDisplaceCenter; cent_x = cent_x + param.param_device.GridSpaceCenter)
210
			{
211
				for (int cent_y = -param.param_device.maxDisplaceCenter; cent_y <= param.param_device.maxDisplaceCenter; cent_y = cent_y + param.param_device.GridSpaceCenter)
212
				{
213
					compareRefMap_kernel<<<divup(maxRef, CUDA_THREAD_COUNT), CUDA_THREAD_COUNT, 0, cudaStream>>> (iOrient, iConv, pConvMap_device[iConv & 1], pProb_device, param.param_device, *pRefMap_device_Mod, cent_x, cent_y, maxRef);
214 215
				}
			}
216
		}
217
		else if (GPUAlgo == 0) //All shifts in one kernel
218
		{
219
			compareRefMapShifted_kernel<<<divup(maxRef, CUDA_THREAD_COUNT), CUDA_THREAD_COUNT, 0, cudaStream>>> (iOrient, iConv, pConvMap_device[iConv & 1], pProb_device, param.param_device, *pRefMap_device_Mod, maxRef);
220
		}
221
		else
222
		{
223 224
			cout << "Invalid GPU Algorithm selected\n";
			exit(1);
225
		}
226
	}
227 228
	if (GPUWorkload < 100)
	{
229
		bioem::compareRefMaps(iOrient, iConv, conv_map, localmultFFT, sumC, sumsquareC, maxRef);
230
	}
231 232 233
	if (GPUAsync)
	{
		checkCudaErrors(cudaEventRecord(cudaEvent[iConv & 1], cudaStream));
234
	}
235 236 237 238 239 240 241 242 243 244
	else
	{
		checkCudaErrors(cudaStreamSynchronize(cudaStream));
	}
	return(0);
}

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

246 247
	if (FFTAlgo) GPUAlgo = 2;

248 249 250 251 252
	gpumap = new bioem_RefMap;
	memcpy(gpumap, &RefMap, sizeof(bioem_RefMap));
	if (FFTAlgo == 0)
	{
		checkCudaErrors(cudaMalloc(&maps, sizeof(myfloat_t) * RefMap.ntotRefMap * RefMap.refMapSize));
253 254 255 256 257 258 259 260 261 262 263 264 265

		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));
		}
266 267 268 269 270 271 272 273 274
	}
	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;

275
	checkCudaErrors(cudaStreamCreate(&cudaStream));
276 277 278
	pProb_device = pProb;
	checkCudaErrors(cudaMalloc(&pProb_device.ptr, pProb_device.get_size(RefMap.ntotRefMap, param.nTotGridAngles)));
	pProb_device.set_pointers();
279
	for (int i = 0; i < 2; i++)
280 281
	{
		checkCudaErrors(cudaEventCreate(&cudaEvent[i]));
282
		checkCudaErrors(cudaMalloc(&pConvMap_device[i], sizeof(myfloat_t) * RefMap.refMapSize));
283
	}
284

285 286
	if (FFTAlgo)
	{
287 288
		checkCudaErrors(cudaMalloc(&pRefMapsFFT, RefMap.ntotRefMap * param.FFTMapSize * sizeof(mycomplex_t)));
		checkCudaErrors(cudaMalloc(&pFFTtmp2, CUDA_FFTS_AT_ONCE * param.FFTMapSize * sizeof(mycomplex_t)));
289
		checkCudaErrors(cudaMalloc(&pFFTtmp, CUDA_FFTS_AT_ONCE * param.param_device.NumberPixels * param.param_device.NumberPixels * sizeof(myfloat_t)));
290 291
		checkCudaErrors(cudaMalloc(&pConvMapFFT, param.FFTMapSize * sizeof(mycomplex_t) * 2));
		checkCudaErrors(cudaMemcpy(pRefMapsFFT, RefMap.RefMapsFFT, RefMap.ntotRefMap * param.FFTMapSize * sizeof(mycomplex_t), cudaMemcpyHostToDevice));
292 293
	}

294 295 296 297 298 299 300
	deviceInitialized = 1;
	return(0);
}

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

302
	cudaStreamDestroy(cudaStream);
303
	cudaFree(pProb_device.ptr);
304 305
	cudaFree(sum);
	cudaFree(sumsquare);
306
	for (int i = 0; i < 2; i++)
307 308
	{
		cudaEventDestroy(cudaEvent[i]);
309
		cudaFree(pConvMap_device[i]);
310
	}
311 312 313 314
	if (FFTAlgo)
	{
		cudaFree(pRefMapsFFT);
		cudaFree(pConvMapFFT);
315
		//cudaFree(pFFTtmp);
316 317
		cudaFree(pFFTtmp2);
	}
318 319 320 321 322 323 324 325 326
	else
	{
		cudaFree(maps);
	}
	if (GPUAlgo == 0 || GPUAlgo == 1)
	{
		cudaFree(pRefMap_device_Mod);
	}
	delete gpumap;
327
	cudaThreadExit();
328

329 330 331 332 333 334
	deviceInitialized = 0;
	return(0);
}

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

337
	cudaMemcpy(pProb_device.ptr, pProb.ptr, pProb.get_size(RefMap.ntotRefMap, param.nTotGridAngles), cudaMemcpyHostToDevice);
338 339 340

	if (FFTAlgo)
	{
341
		for (int i = 0; i < 2; i++)
342 343
		{
			int n[2] = {param.param_device.NumberPixels, param.param_device.NumberPixels};
344
			if (cufftPlanMany(&plan[i], 2, n, NULL, 1, param.FFTMapSize, NULL, 1, 0, MY_CUFFT_C2R, i ? (maxRef % CUDA_FFTS_AT_ONCE) : CUDA_FFTS_AT_ONCE) != CUFFT_SUCCESS)
345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360
			{
				cout << "Error planning CUFFT\n";
				exit(1);
			}
			if (cufftSetCompatibilityMode(plan[i], CUFFT_COMPATIBILITY_NATIVE) != CUFFT_SUCCESS)
			{
				cout << "Error planning CUFFT compatibility\n";
				exit(1);
			}
			if (cufftSetStream(plan[i], cudaStream) != CUFFT_SUCCESS)
			{
				cout << "Error setting CUFFT stream\n";
				exit(1);
			}
		}
	}
361 362 363 364 365 366
	return(0);
}

int bioem_cuda::deviceFinishRun()
{
	if (GPUAsync) cudaStreamSynchronize(cudaStream);
367
	cudaMemcpy(pProb.ptr, pProb_device.ptr, pProb.get_size(RefMap.ntotRefMap, param.nTotGridAngles), cudaMemcpyDeviceToHost);
368

369 370
	if (FFTAlgo)
	{
371
		for (int i = 0; i < 2; i++) cufftDestroy(plan[i]);
372 373
	}

374 375 376 377 378 379 380
	return(0);
}

bioem* bioem_cuda_create()
{
	return new bioem_cuda;
}