bioem_cuda.cu 27.2 KB
Newer Older
Pilar Cossio's avatar
License  
Pilar Cossio committed
1
/* ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 3 4 5
   < BioEM software for Bayesian inference of Electron Microscopy images>
   Copyright (C) 2016 Pilar Cossio, David Rohr, Fabio Baruffa, Markus Rampp, 
        Volker Lindenstruth and Gerhard Hummer.
   Max Planck Institute of Biophysics, Frankfurt, Germany.
6

7
   See license statement for terms of distribution.
Pilar Cossio's avatar
License  
Pilar Cossio committed
8 9 10

   ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++*/

11 12 13 14 15 16 17 18 19 20 21
#define BIOEM_GPUCODE

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

#include <iostream>
using namespace std;

#include "bioem_cuda_internal.h"
#include "bioem_algorithm.h"
Pilar Cossio's avatar
Pilar Cossio committed
22
//#include "helper_cuda.h"
23

24 25 26 27 28 29 30
#define checkCudaErrors(error) \
{ \
	if ((error) != cudaSuccess) \
	{ \
		printf("CUDA Error %d / %s (%s: %d)\n", error, cudaGetErrorString(error), __FILE__, __LINE__); \
		exit(1); \
	} \
31 32
}

David Rohr's avatar
David Rohr committed
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
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";
}

Luka Stanisic's avatar
Luka Stanisic committed
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
/* Handing CUDA Driver errors */

#define cuErrorCheck(call) \
  do { \
    CUresult __error__; \
    if ((__error__ = (call)) != CUDA_SUCCESS) { \
      printf("CUDA Driver Error %d / %s (%s %d)\n", __error__, cuGetError(__error__),__FILE__, __LINE__); \
      return __error__; \
    } \
  } while (false)

static const char * cuGetError(CUresult result) {
  switch (result) {
    case CUDA_SUCCESS:                              return "No errors";
    case CUDA_ERROR_INVALID_VALUE:                  return "Invalid value";
    case CUDA_ERROR_OUT_OF_MEMORY:                  return "Out of memory";
    case CUDA_ERROR_NOT_INITIALIZED:                return "Driver not initialized";
    case CUDA_ERROR_DEINITIALIZED:                  return "Driver deinitialized";
    case CUDA_ERROR_PROFILER_DISABLED:              return "Profiler disabled";
    case CUDA_ERROR_PROFILER_NOT_INITIALIZED:       return "Profiler not initialized";
    case CUDA_ERROR_PROFILER_ALREADY_STARTED:       return "Profiler already started";
    case CUDA_ERROR_PROFILER_ALREADY_STOPPED:       return "Profiler already stopped";
    case CUDA_ERROR_NO_DEVICE:                      return "No CUDA-capable device available";
    case CUDA_ERROR_INVALID_DEVICE:                 return "Invalid device";
    case CUDA_ERROR_INVALID_IMAGE:                  return "Invalid kernel image";
    case CUDA_ERROR_INVALID_CONTEXT:                return "Invalid context";
    case CUDA_ERROR_CONTEXT_ALREADY_CURRENT:        return "Context already current";
    case CUDA_ERROR_MAP_FAILED:                     return "Map failed";
    case CUDA_ERROR_UNMAP_FAILED:                   return "Unmap failed";
    case CUDA_ERROR_ARRAY_IS_MAPPED:                return "Array is mapped";
    case CUDA_ERROR_ALREADY_MAPPED:                 return "Already mapped";
    case CUDA_ERROR_NO_BINARY_FOR_GPU:              return "No binary for GPU";
    case CUDA_ERROR_ALREADY_ACQUIRED:               return "Already acquired";
    case CUDA_ERROR_NOT_MAPPED:                     return "Not mapped";
    case CUDA_ERROR_NOT_MAPPED_AS_ARRAY:            return "Not mapped as array";
    case CUDA_ERROR_NOT_MAPPED_AS_POINTER:          return "Not mapped as pointer";
    case CUDA_ERROR_ECC_UNCORRECTABLE:              return "Uncorrectable ECC error";
    case CUDA_ERROR_UNSUPPORTED_LIMIT:              return "Unsupported CUlimit";
    case CUDA_ERROR_CONTEXT_ALREADY_IN_USE:         return "Context already in use";
    case CUDA_ERROR_INVALID_SOURCE:                 return "Invalid source";
    case CUDA_ERROR_FILE_NOT_FOUND:                 return "File not found";
    case CUDA_ERROR_SHARED_OBJECT_SYMBOL_NOT_FOUND: return "Shared object symbol not found";
    case CUDA_ERROR_SHARED_OBJECT_INIT_FAILED:      return "Shared object initialization failed";
    case CUDA_ERROR_OPERATING_SYSTEM:               return "Operating System call failed";
    case CUDA_ERROR_INVALID_HANDLE:                 return "Invalid handle";
    case CUDA_ERROR_NOT_FOUND:                      return "Not found";
    case CUDA_ERROR_NOT_READY:                      return "CUDA not ready";
    case CUDA_ERROR_LAUNCH_FAILED:                  return "Launch failed";
    case CUDA_ERROR_LAUNCH_OUT_OF_RESOURCES:        return "Launch exceeded resources";
    case CUDA_ERROR_LAUNCH_TIMEOUT:                 return "Launch exceeded timeout";
    case CUDA_ERROR_LAUNCH_INCOMPATIBLE_TEXTURING:  return "Launch with incompatible texturing";
    case CUDA_ERROR_PEER_ACCESS_ALREADY_ENABLED:    return "Peer access already enabled";
    case CUDA_ERROR_PEER_ACCESS_NOT_ENABLED:        return "Peer access not enabled";
    case CUDA_ERROR_PRIMARY_CONTEXT_ACTIVE:         return "Primary context active";
    case CUDA_ERROR_CONTEXT_IS_DESTROYED:           return "Context is destroyed";
    case CUDA_ERROR_ASSERT:                         return "Device assert failed";
    case CUDA_ERROR_TOO_MANY_PEERS:                 return "Too many peers";
    case CUDA_ERROR_HOST_MEMORY_ALREADY_REGISTERED: return "Host memory already registered";
    case CUDA_ERROR_HOST_MEMORY_NOT_REGISTERED:     return "Host memory not registered";
    case CUDA_ERROR_UNKNOWN:                        return "Unknown error";
    default:                                        return "Unknown error code";
  }
}

134 135 136 137
bioem_cuda::bioem_cuda()
{
	deviceInitialized = 0;
	GPUAlgo = getenv("GPUALGO") == NULL ? 2 : atoi(getenv("GPUALGO"));
138 139
	GPUAsync = getenv("GPUASYNC") == NULL ? 1 : atoi(getenv("GPUASYNC"));
	GPUWorkload = getenv("GPUWORKLOAD") == NULL ? 100 : atoi(getenv("GPUWORKLOAD"));
Luka Stanisic's avatar
Luka Stanisic committed
140
	if (GPUWorkload == -1) GPUWorkload = 100;
141
	GPUDualStream = getenv("GPUDUALSTREAM") == NULL ? 1 : atoi(getenv("GPUDUALSTREAM"));
142 143 144 145 146 147 148
}

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

149 150 151
__global__ void compareRefMap_kernel(const int iOrient, const int iConv,  const myfloat_t amp, const myfloat_t pha, const myfloat_t env, const myfloat_t sumC,
                                                const myfloat_t sumsquareC, 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)
152 153
{
	const int iRefMap = myBlockIdxX * myBlockDimX + myThreadIdxX;
154
	if (iRefMap < maxRef)
155
	{
156
		compareRefMap<0>(iRefMap, iOrient, iConv, amp, pha, env, sumC, sumsquareC, pMap, pProb, param, RefMap, cent_x, cent_y);
157 158 159
	}
}

Pilar Cossio's avatar
Pilar Cossio committed
160
__global__ void compareRefMapShifted_kernel(const int iOrient, const int iConv, const myfloat_t amp, const myfloat_t pha, const myfloat_t env, const myfloat_t sumC, const myfloat_t sumsquareC, const myfloat_t* pMap, bioem_Probability pProb, const bioem_param_device param, const bioem_RefMap_Mod RefMap, const int maxRef)
161 162
{
	const int iRefMap = myBlockIdxX * myBlockDimX + myThreadIdxX;
163
	if (iRefMap < maxRef)
164
	{
165
		compareRefMapShifted<1>(iRefMap, iOrient, iConv, amp, pha, env, sumC, sumsquareC, pMap, pProb, param, RefMap);
166 167 168
	}
}

169 170 171 172 173 174
__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;
175
	for (int i = myid; i < mysize; i += mygrid) myptr[i] = 0;
176 177
}

Pilar Cossio's avatar
Pilar Cossio committed
178
__global__ void compareRefMapLoopShifts_kernel(const int iOrient, const int iConv, const myfloat_t amp, const myfloat_t pha, const myfloat_t env, const myfloat_t sumC, const myfloat_t sumsquareC, 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)
179 180 181 182 183 184 185 186 187
{
	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;
188

189
	const bool threadActive = myShiftIdx < nShifts && myShiftIdy < nShifts && iRefMap < maxRef;
190

Pilar Cossio's avatar
Pilar Cossio committed
191
	compareRefMap<2>(iRefMap, iOrient, iConv, amp, pha, env, sumC, sumsquareC, pMap, pProb, param, RefMap, cent_x, cent_y, myShift, nShifts * nShifts, myRef, threadActive);
192 193
}

194
__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)
195 196
{
	if (myBlockIdxX >= NumberMaps) return;
197
	const mycuComplex_t* myin = (mycuComplex_t*) &refmap[(myBlockIdxX + Offset) * MapSize];
198
	const mycuComplex_t* myconv = (mycuComplex_t*) convmap;
199
	mycuComplex_t* myout = &out[myBlockIdxX * MapSize];
200
	for(int i = myThreadIdxX; i < NumberPixelsTotal; i += myBlockDimX)
201
	{
202 203 204 205 206 207 208
		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;
209 210 211
	}
}

212
__global__ void cuDoRefMapsFFT(const int iOrient, const int iConv, const myfloat_t amp, const myfloat_t pha, const myfloat_t env, 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)
213
{
214
	if (myBlockIdxX * myBlockDimX + myThreadIdxX >= maxRef) return;
215 216
	const int iRefMap = myBlockIdxX * myBlockDimX + myThreadIdxX + Offset;
	const myfloat_t* mylCC = &lCC[(myBlockIdxX * myBlockDimX + myThreadIdxX) * param.NumberPixels * param.NumberPixels];
217
	doRefMapFFT(iRefMap, iOrient, iConv, amp, pha, env, mylCC, sumC, sumsquareC, pProb, param, RefMap);
218 219
}

220 221 222 223 224 225 226 227 228 229 230 231 232
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

233
int bioem_cuda::compareRefMaps(int iOrient, int iConv, myfloat_t amp, myfloat_t pha, myfloat_t env, const myfloat_t* conv_map, mycomplex_t* localmultFFT, myfloat_t sumC, myfloat_t sumsquareC, const int startMap)
234
{
235 236 237 238 239
	if (startMap)
	{
		cout << "Error startMap not implemented for GPU Code\n";
		exit(1);
	}
Luka Stanisic's avatar
Luka Stanisic committed
240 241 242 243 244 245 246
#ifdef DEBUG_GPU
	float time;
	cudaEvent_t start, stop;
	checkCudaErrors(cudaEventCreate(&start));
	checkCudaErrors(cudaEventCreate(&stop));
	checkCudaErrors(cudaEventRecord(start, 0));
#endif
247 248 249 250
	if (GPUAsync)
	{
		checkCudaErrors(cudaEventSynchronize(cudaEvent[iConv & 1]));
	}
Luka Stanisic's avatar
Luka Stanisic committed
251 252 253 254 255 256 257
#ifdef DEBUG_GPU
	checkCudaErrors(cudaEventRecord(stop, 0));
	checkCudaErrors(cudaEventSynchronize(stop));
	checkCudaErrors(cudaEventElapsedTime(&time, start, stop));
	printf("\t\t\tGPU: time to synch projections %1.6f sec\n", time/1000);
	checkCudaErrors(cudaEventRecord(start, 0));
#endif
258
	if (FFTAlgo)
259
	{
260
		memcpy(&pConvMapFFT_Host[(iConv & 1) * param.FFTMapSize], localmultFFT, param.FFTMapSize * sizeof(mycomplex_t));
261
		checkCudaErrors(cudaMemcpyAsync(&pConvMapFFT[(iConv & 1) * param.FFTMapSize], &pConvMapFFT_Host[(iConv & 1) * param.FFTMapSize], param.FFTMapSize * sizeof(mycomplex_t), cudaMemcpyHostToDevice, cudaStream[GPUAsync ? 2 : 0]));
Luka Stanisic's avatar
Luka Stanisic committed
262 263 264 265 266 267 268
#ifdef DEBUG_GPU
		checkCudaErrors(cudaEventRecord(stop, 0));
		checkCudaErrors(cudaEventSynchronize(stop));
		checkCudaErrors(cudaEventElapsedTime(&time, start, stop));
		printf("\t\t\tGPU: time for memcpy %1.6f sec\n", time/1000);
		checkCudaErrors(cudaEventRecord(start, 0));
#endif
269 270 271 272 273
		if (GPUAsync)
		{
			checkCudaErrors(cudaEventRecord(cudaEvent[2], cudaStream[2]));
			checkCudaErrors(cudaStreamWaitEvent(cudaStream[0], cudaEvent[2], 0));
		}
274
		if (GPUDualStream)
275
		{
276 277 278 279 280 281
			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;
282
			const int num = min(CUDA_FFTS_AT_ONCE, maxRef - i);
283 284
			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
285
			if (err != CUFFT_SUCCESS)
286
			{
David Rohr's avatar
David Rohr committed
287
				cout << "Error running CUFFT " << cufftGetErrorStrung(err) << "\n";
288 289
				exit(1);
			}
290
			cuDoRefMapsFFT<<<divup(num, CUDA_THREAD_COUNT), CUDA_THREAD_COUNT, 0, cudaStream[j & 1]>>>(iOrient, iConv,  amp, pha, env, pFFTtmp[j & 1], sumC, sumsquareC, pProb_device, param.param_device, *gpumap, num, i);
291
		}
Luka Stanisic's avatar
Luka Stanisic committed
292
		checkCudaErrors(cudaPeekAtLastError());
293 294 295 296 297
		if (GPUDualStream)
		{
			checkCudaErrors(cudaEventRecord(cudaFFTEvent[1], cudaStream[1]));
			checkCudaErrors(cudaStreamWaitEvent(cudaStream[0], cudaFFTEvent[1], 0));
		}
298 299 300
	}
	else
	{
301
		checkCudaErrors(cudaMemcpyAsync(pConvMap_device[iConv & 1], conv_map, sizeof(myfloat_t) * RefMap.refMapSize, cudaMemcpyHostToDevice, cudaStream[0]));
Luka Stanisic's avatar
Luka Stanisic committed
302 303 304 305 306 307 308
#ifdef DEBUG_GPU
		checkCudaErrors(cudaEventRecord(stop, 0));
		checkCudaErrors(cudaEventSynchronize(stop));
		checkCudaErrors(cudaEventElapsedTime(&time, start, stop));
		printf("\t\t\tGPU: time for memcpy %1.6f sec\n", time/1000);
		checkCudaErrors(cudaEventRecord(start, 0) );
#endif
309
		if (GPUAlgo == 2) //Loop over shifts
310
		{
311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329
			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;
330
			for (size_t i = 0; i < totalBlocks; i += nBlocks)
331
			{
Pilar Cossio's avatar
Pilar Cossio committed
332
				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, amp, pha, env, sumC, sumsquareC, pConvMap_device[iConv & 1], pProb_device, param.param_device, *gpumap, i, nShifts, nShiftBits, maxRef);
333
			}
334
		}
335
		else if (GPUAlgo == 1) //Split shifts in multiple kernels
336
		{
337
			for (int cent_x = -param.param_device.maxDisplaceCenter; cent_x <= param.param_device.maxDisplaceCenter; cent_x = cent_x + param.param_device.GridSpaceCenter)
338
			{
339
				for (int cent_y = -param.param_device.maxDisplaceCenter; cent_y <= param.param_device.maxDisplaceCenter; cent_y = cent_y + param.param_device.GridSpaceCenter)
340
				{
Pilar Cossio's avatar
Pilar Cossio committed
341
					compareRefMap_kernel<<<divup(maxRef, CUDA_THREAD_COUNT), CUDA_THREAD_COUNT, 0, cudaStream[0]>>> (iOrient, iConv, amp, pha, env, sumC, sumsquareC, pConvMap_device[iConv & 1], pProb_device, param.param_device, *pRefMap_device_Mod, cent_x, cent_y, maxRef);
342 343
				}
			}
344
		}
345
		else if (GPUAlgo == 0) //All shifts in one kernel
346
		{
347
			compareRefMapShifted_kernel<<<divup(maxRef, CUDA_THREAD_COUNT), CUDA_THREAD_COUNT, 0, cudaStream[0]>>> (iOrient, iConv, amp, pha, env, sumC, sumsquareC, pConvMap_device[iConv & 1], pProb_device, param.param_device, *pRefMap_device_Mod, maxRef);
348
		}
349
		else
350
		{
351 352
			cout << "Invalid GPU Algorithm selected\n";
			exit(1);
353
		}
354
	}
Luka Stanisic's avatar
Luka Stanisic committed
355 356 357 358 359 360 361
#ifdef DEBUG_GPU
	checkCudaErrors(cudaEventRecord(stop, 0));
	checkCudaErrors(cudaEventSynchronize(stop));
	checkCudaErrors(cudaEventElapsedTime(&time, start, stop));
	printf("\t\t\tGPU: time to run CUDA %1.6f sec\n", time/1000);
	checkCudaErrors(cudaEventRecord(start, 0));
#endif
362 363
	if (GPUWorkload < 100)
	{
364
		bioem::compareRefMaps(iOrient, iConv, amp, pha, env, conv_map, localmultFFT, sumC, sumsquareC, maxRef);
365
	}
Luka Stanisic's avatar
Luka Stanisic committed
366 367 368 369 370 371
#ifdef DEBUG_GPU
	checkCudaErrors(cudaEventRecord(stop, 0));
	checkCudaErrors(cudaEventSynchronize(stop));
	checkCudaErrors(cudaEventElapsedTime(&time, start, stop));
	printf("\t\t\tGPU: time to run OMP %1.6f sec\n", time/1000);
#endif
372 373
	if (GPUAsync)
	{
374
		checkCudaErrors(cudaEventRecord(cudaEvent[iConv & 1], cudaStream[0]));
375
	}
376 377
	else
	{
378
		checkCudaErrors(cudaStreamSynchronize(cudaStream[0]));
379 380 381 382
	}
	return(0);
}

David Rohr's avatar
David Rohr committed
383 384 385
int bioem_cuda::selectCudaDevice()
{
	int count;
386
	int bestDevice = 0;
David Rohr's avatar
David Rohr committed
387
	cudaDeviceProp deviceProp;
388 389 390 391 392

	/* Initializing CUDA driver API */
	cuErrorCheck(cuInit(0));

	/* Get number of available CUDA devices */
David Rohr's avatar
David Rohr committed
393 394 395 396 397 398
	checkCudaErrors(cudaGetDeviceCount(&count));
	if (count == 0)
	{
		printf("No CUDA device detected\n");
		return(1);
	}
399 400 401 402

	/* Find the best GPU */
	long long int bestDeviceSpeed = -1, deviceSpeed = -1;
	for (int i = 0; i < count; i++)
David Rohr's avatar
David Rohr committed
403
	{
404 405
		cudaGetDeviceProperties(&deviceProp, i);
		deviceSpeed = (long long int) deviceProp.multiProcessorCount * (long long int) deviceProp.clockRate * (long long int) deviceProp.warpSize;
David Rohr's avatar
David Rohr committed
406 407 408 409 410 411
		if (deviceSpeed > bestDeviceSpeed)
		{
			bestDevice = i;
			bestDeviceSpeed = deviceSpeed;
		}
	}
412 413

	/* Get user-specified GPU choice */
414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430
	if (getenv("GPUDEVICE"))
	{
		int device = atoi(getenv("GPUDEVICE"));
		if (device > count)
		{
			printf("Invalid CUDA device specified, max device number is %d\n", count);
			exit(1);
		}
#ifdef WITH_MPI
		if (device == -1)
		{
			device = mpi_rank % count;
		}
#endif
		if (device < 0)
		{
			printf("Negative CUDA device specified: %d, invalid!\n", device);
431
			exit(1);
432 433 434
		}
		bestDevice = device;
	}
David Rohr's avatar
David Rohr committed
435

436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452
	/* Set CUDA processes to appropriate devices */
	cudaGetDeviceProperties(&deviceProp, bestDevice);
	if (deviceProp.computeMode == 0)
	{
		checkCudaErrors(cudaSetDevice(bestDevice));
	}
	else
	{
		if (DebugOutput >= 1)
		{
			printf("CUDA device %d is not set in DEFAULT mode, make sure that CUDA processes are pinned as planned!\n", bestDevice);
			printf("Pinning process %d to CUDA device %d\n", mpi_rank, bestDevice);
		}
		checkCudaErrors(cudaSetDevice(bestDevice));
		/* This synchronization is needed in order to detect bogus silent errors from cudaSetDevice call */
		checkCudaErrors(cudaDeviceSynchronize());
	}
David Rohr's avatar
David Rohr committed
453

454
	/* Debugging information about CUDA devices used by the current process */
David Rohr's avatar
David Rohr committed
455
	if (DebugOutput >= 3)
David Rohr's avatar
David Rohr committed
456
	{
David Rohr's avatar
David Rohr committed
457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472
		printf("Using CUDA Device %s with Properties:\n", deviceProp.name);
		printf("totalGlobalMem = %lld\n", (unsigned long long int) deviceProp.totalGlobalMem);
		printf("sharedMemPerBlock = %lld\n", (unsigned long long int) deviceProp.sharedMemPerBlock);
		printf("regsPerBlock = %d\n", deviceProp.regsPerBlock);
		printf("warpSize = %d\n", deviceProp.warpSize);
		printf("memPitch = %lld\n", (unsigned long long int) deviceProp.memPitch);
		printf("maxThreadsPerBlock = %d\n", deviceProp.maxThreadsPerBlock);
		printf("maxThreadsDim = %d %d %d\n", deviceProp.maxThreadsDim[0], deviceProp.maxThreadsDim[1], deviceProp.maxThreadsDim[2]);
		printf("maxGridSize = %d %d %d\n", deviceProp.maxGridSize[0], deviceProp.maxGridSize[1], deviceProp.maxGridSize[2]);
		printf("totalConstMem = %lld\n", (unsigned long long int) deviceProp.totalConstMem);
		printf("major = %d\n", deviceProp.major);
		printf("minor = %d\n", deviceProp.minor);
		printf("clockRate = %d\n", deviceProp.clockRate);
		printf("memoryClockRate = %d\n", deviceProp.memoryClockRate);
		printf("multiProcessorCount = %d\n", deviceProp.multiProcessorCount);
		printf("textureAlignment = %lld\n", (unsigned long long int) deviceProp.textureAlignment);
473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492
		printf("computeMode = %d\n", deviceProp.computeMode);
#if CUDA_VERSION > 3010
		size_t free, total;
#else
		unsigned int free, total;
#endif
		if (deviceProp.computeMode == 0)
		{
			CUdevice tmpDevice;
			cuErrorCheck(cuDeviceGet(&tmpDevice, bestDevice));
			CUcontext tmpContext;
			cuErrorCheck(cuCtxCreate(&tmpContext, 0, tmpDevice));
			cuErrorCheck(cuMemGetInfo(&free, &total));
			cuErrorCheck(cuCtxDestroy(tmpContext));
		}
		else
		{
			cuErrorCheck(cuMemGetInfo(&free, &total));
		}
		printf("free memory = %lld; total memory = %lld\n", free, total);
David Rohr's avatar
David Rohr committed
493
	}
494

David Rohr's avatar
David Rohr committed
495 496
	if (DebugOutput >= 1)
	{
David Rohr's avatar
David Rohr committed
497
		printf("BioEM for CUDA initialized (MPI Rank %d), %d GPUs found, using GPU %d\n", mpi_rank, count, bestDevice);
David Rohr's avatar
David Rohr committed
498
	}
499

David Rohr's avatar
David Rohr committed
500 501 502
	return(0);
}

503 504 505
int bioem_cuda::deviceInit()
{
	deviceExit();
David Rohr's avatar
David Rohr committed
506
	
507
	selectCudaDevice();
508

509 510
	if (FFTAlgo) GPUAlgo = 2;

511 512 513 514 515
	gpumap = new bioem_RefMap;
	memcpy(gpumap, &RefMap, sizeof(bioem_RefMap));
	if (FFTAlgo == 0)
	{
		checkCudaErrors(cudaMalloc(&maps, sizeof(myfloat_t) * RefMap.ntotRefMap * RefMap.refMapSize));
516 517 518 519 520 521 522 523 524 525 526 527 528

		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));
		}
529 530 531 532 533 534 535 536 537
	}
	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;

538
	checkCudaErrors(cudaMalloc(&pProb_memory, pProb_device.get_size(RefMap.ntotRefMap, param.nTotGridAngles, param.nTotCC, param.param_device.writeAngles, param.param_device.writeCC)));
539

540
	for (int i = 0; i < 2; i++)
541
	{
542
		checkCudaErrors(cudaStreamCreate(&cudaStream[i]));
543
		checkCudaErrors(cudaEventCreate(&cudaEvent[i]));
544
		checkCudaErrors(cudaEventCreate(&cudaFFTEvent[i]));
545
		checkCudaErrors(cudaMalloc(&pConvMap_device[i], sizeof(myfloat_t) * RefMap.refMapSize));
546
	}
547 548 549 550 551
	if (GPUAsync)
	{
		checkCudaErrors(cudaStreamCreate(&cudaStream[2]));
		checkCudaErrors(cudaEventCreate(&cudaEvent[2]));
	}
552

553 554
	if (FFTAlgo)
	{
555
		checkCudaErrors(cudaMalloc(&pRefMapsFFT, RefMap.ntotRefMap * param.FFTMapSize * sizeof(mycomplex_t)));
556 557 558 559
		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;
560
		checkCudaErrors(cudaMalloc(&pConvMapFFT, param.FFTMapSize * sizeof(mycomplex_t) * 2));
561
		checkCudaErrors(cudaHostAlloc(&pConvMapFFT_Host, param.FFTMapSize * sizeof(mycomplex_t) * 2, 0));
562
		checkCudaErrors(cudaMemcpy(pRefMapsFFT, RefMap.RefMapsFFT, RefMap.ntotRefMap * param.FFTMapSize * sizeof(mycomplex_t), cudaMemcpyHostToDevice));
563 564
	}

565 566 567 568 569 570 571
	deviceInitialized = 1;
	return(0);
}

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

573

David Rohr's avatar
David Rohr committed
574
	cudaFree(pProb_memory);
575 576
	cudaFree(sum);
	cudaFree(sumsquare);
577
	for (int i = 0; i < 2; i++)
578
	{
579
		cudaStreamDestroy(cudaStream[i]);
580
		cudaEventDestroy(cudaEvent[i]);
581
		cudaEventDestroy(cudaFFTEvent[i]);
582
		cudaFree(pConvMap_device[i]);
583
	}
584 585 586 587
	if (FFTAlgo)
	{
		cudaFree(pRefMapsFFT);
		cudaFree(pConvMapFFT);
588
		cudaFreeHost(pConvMapFFT_Host);
589 590
		cudaFree(pFFTtmp[0]);
		cudaFree(pFFTtmp2[0]);
591
	}
592 593 594 595 596 597 598 599
	else
	{
		cudaFree(maps);
	}
	if (GPUAlgo == 0 || GPUAlgo == 1)
	{
		cudaFree(pRefMap_device_Mod);
	}
600 601 602 603 604 605
	if (GPUAsync)
	{
		cudaStreamDestroy(cudaStream[2]);
		cudaEventDestroy(cudaEvent[2]);
	}

606
	delete gpumap;
607
	cudaThreadExit();
608

609 610 611 612 613 614
	deviceInitialized = 0;
	return(0);
}

int bioem_cuda::deviceStartRun()
{
David Rohr's avatar
David Rohr committed
615 616 617 618 619 620 621
	if (GPUWorkload >= 100)
	{
		maxRef = RefMap.ntotRefMap;
		pProb_host = &pProb;
	}
	else
	{
622
		maxRef = RefMap.ntotRefMap == 1 ? (size_t) RefMap.ntotRefMap : (size_t) RefMap.ntotRefMap * (size_t) GPUWorkload / 100;
David Rohr's avatar
David Rohr committed
623
		pProb_host = new bioem_Probability;
624
		pProb_host->init(maxRef, param.nTotGridAngles, param.nTotCC, *this);
David Rohr's avatar
David Rohr committed
625 626
		pProb_host->copyFrom(&pProb, *this);
	}
627

David Rohr's avatar
David Rohr committed
628 629 630
	pProb_device = *pProb_host;
	pProb_device.ptr = pProb_memory;
	pProb_device.set_pointers();
631
	checkCudaErrors(cudaMemcpyAsync(pProb_device.ptr, pProb_host->ptr, pProb_host->get_size(maxRef, param.nTotGridAngles, param.nTotCC, param.param_device.writeAngles, param.param_device.writeCC), cudaMemcpyHostToDevice, cudaStream[0]));
632 633 634

	if (FFTAlgo)
	{
635
		for (int j = 0;j < 2;j++)
636
		{
637
			for (int i = 0; i < 2; i++)
638
			{
639
				if (i && maxRef % CUDA_FFTS_AT_ONCE == 0) continue;
640 641 642 643 644 645
				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);
				}
646
			        if (cufftSetCompatibilityMode(plan[i][j], CUFFT_COMPATIBILITY_FFTW_PADDING) != CUFFT_SUCCESS)
647 648 649 650 651 652 653 654 655
				{
					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);
				}
656
			}
657
			if (!GPUDualStream) break;
658 659
		}
	}
660 661 662 663 664
	return(0);
}

int bioem_cuda::deviceFinishRun()
{
665
	if (GPUAsync) cudaStreamSynchronize(cudaStream[0]);
666
	checkCudaErrors(cudaMemcpyAsync(pProb_host->ptr, pProb_device.ptr, pProb_host->get_size(maxRef, param.nTotGridAngles, param.nTotCC, param.param_device.writeAngles, param.param_device.writeCC), cudaMemcpyDeviceToHost, cudaStream[0]));
667

668 669
	if (FFTAlgo)
	{
670 671
		for (int j = 0;j < 2;j++)
		{
672 673 674 675 676
			for (int i = 0; i < 2; i++)
			{
				if (i && maxRef % CUDA_FFTS_AT_ONCE == 0) continue;
				cufftDestroy(plan[i][j]);
			}
677 678
			if (!GPUDualStream) break;
		}
679
	}
David Rohr's avatar
David Rohr committed
680 681 682 683
	cudaThreadSynchronize();
	if (GPUWorkload < 100)
	{
		pProb.copyFrom(pProb_host, *this);
684
		free_device_host(pProb_host->ptr);
David Rohr's avatar
David Rohr committed
685 686
		delete[] pProb_host;
	}
687

688 689 690
	return(0);
}

691 692 693 694 695 696 697 698 699 700 701 702
void* bioem_cuda::malloc_device_host(size_t size)
{
	void* ptr;
	checkCudaErrors(cudaHostAlloc(&ptr, size, 0));
	return(ptr);
}

void bioem_cuda::free_device_host(void* ptr)
{
	cudaFreeHost(ptr);
}

Luka Stanisic's avatar
Luka Stanisic committed
703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719
void bioem_cuda::rebalance(int workload)
{
	if ((workload < 0) || (workload > 100) || (workload == GPUWorkload)) return;

	deviceFinishRun();

	if (DebugOutput >= 2)
	{
	  printf("\t\tSetting GPU workload to %d%% (rank %d)\n", workload, mpi_rank);
	}

	GPUWorkload = workload;
	maxRef = (size_t) RefMap.ntotRefMap * (size_t) GPUWorkload / 100;

	deviceStartRun();
}

720 721
bioem* bioem_cuda_create()
{
David Rohr's avatar
David Rohr committed
722 723 724 725 726 727 728 729 730
	int count;
	
	if (cudaGetDeviceCount(&count) != cudaSuccess) count = 0;
	if (count == 0)
	{
		printf("No CUDA device available, using fallback to CPU version\n");
		return new bioem;
	}

731 732
	return new bioem_cuda;
}