bioem_cuda.cu 38.3 KB
Newer Older
Pilar Cossio's avatar
License  
Pilar Cossio committed
1
/* ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2
   < BioEM software for Bayesian inference of Electron Microscopy images>
Luka Stanisic's avatar
Luka Stanisic committed
3 4
   Copyright (C) 2017 Pilar Cossio, David Rohr, Fabio Baruffa, Markus Rampp,
        Luka Stanisic, Volker Lindenstruth and Gerhard Hummer.
5
   Max Planck Institute of Biophysics, Frankfurt, Germany.
Luka Stanisic's avatar
Luka Stanisic committed
6 7 8
   Frankfurt Institute for Advanced Studies, Goethe University Frankfurt,
   Germany.
   Max Planck Computing and Data Facility, Garching, Germany.
9

Luka Stanisic's avatar
Luka Stanisic committed
10
   Released under the GNU Public License, v3.
11
   See license statement for terms of distribution.
Pilar Cossio's avatar
License  
Pilar Cossio committed
12 13 14

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

15 16 17 18 19 20 21 22 23 24
#define BIOEM_GPUCODE

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

#include <iostream>
using namespace std;

#include "bioem_cuda_internal.h"
Pilar Cossio's avatar
Pilar Cossio committed
25
//#include "helper_cuda.h"
26

Luka Stanisic's avatar
Luka Stanisic committed
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
#include "bioem_algorithm.h"

#define checkCudaErrors(error)                                                 \
  {                                                                            \
    if ((error) != cudaSuccess)                                                \
    {                                                                          \
      printf("CUDA Error %d / %s (%s: %d)\n", error,                           \
             cudaGetErrorString(error), __FILE__, __LINE__);                   \
      exit(1);                                                                 \
    }                                                                          \
  }

#ifdef DEBUG_GPU
#define printCudaDebugStart()                                                  \
  float time;                                                                  \
  time = 0.;                                                                   \
  cudaEvent_t start, stop;                                                     \
  checkCudaErrors(cudaEventCreate(&start));                                    \
  checkCudaErrors(cudaEventCreate(&stop));                                     \
  checkCudaErrors(cudaEventRecord(start, 0));
#define printCudaDebug(msg)                                                    \
  checkCudaErrors(cudaEventRecord(stop, 0));                                   \
  checkCudaErrors(cudaEventSynchronize(stop));                                 \
  checkCudaErrors(cudaEventElapsedTime(&time, start, stop));                   \
  printf("\t\t\tGPU: %s %1.6f sec\n", msg, time / 1000);                       \
  checkCudaErrors(cudaEventRecord(start, 0));

#else
#define printCudaDebugStart()
#define printCudaDebug(msg)
#endif
58

David Rohr's avatar
David Rohr committed
59 60
static const char *cufftGetErrorStrung(cufftResult error)
{
Luka Stanisic's avatar
Luka Stanisic committed
61 62 63 64
  switch (error)
  {
    case CUFFT_SUCCESS:
      return "CUFFT_SUCCESS";
David Rohr's avatar
David Rohr committed
65

Luka Stanisic's avatar
Luka Stanisic committed
66 67
    case CUFFT_INVALID_PLAN:
      return "CUFFT_INVALID_PLAN";
David Rohr's avatar
David Rohr committed
68

Luka Stanisic's avatar
Luka Stanisic committed
69 70
    case CUFFT_ALLOC_FAILED:
      return "CUFFT_ALLOC_FAILED";
David Rohr's avatar
David Rohr committed
71

Luka Stanisic's avatar
Luka Stanisic committed
72 73
    case CUFFT_INVALID_TYPE:
      return "CUFFT_INVALID_TYPE";
David Rohr's avatar
David Rohr committed
74

Luka Stanisic's avatar
Luka Stanisic committed
75 76
    case CUFFT_INVALID_VALUE:
      return "CUFFT_INVALID_VALUE";
David Rohr's avatar
David Rohr committed
77

Luka Stanisic's avatar
Luka Stanisic committed
78 79
    case CUFFT_INTERNAL_ERROR:
      return "CUFFT_INTERNAL_ERROR";
David Rohr's avatar
David Rohr committed
80

Luka Stanisic's avatar
Luka Stanisic committed
81 82
    case CUFFT_EXEC_FAILED:
      return "CUFFT_EXEC_FAILED";
David Rohr's avatar
David Rohr committed
83

Luka Stanisic's avatar
Luka Stanisic committed
84 85
    case CUFFT_SETUP_FAILED:
      return "CUFFT_SETUP_FAILED";
David Rohr's avatar
David Rohr committed
86

Luka Stanisic's avatar
Luka Stanisic committed
87 88
    case CUFFT_INVALID_SIZE:
      return "CUFFT_INVALID_SIZE";
David Rohr's avatar
David Rohr committed
89

Luka Stanisic's avatar
Luka Stanisic committed
90 91 92 93
    case CUFFT_UNALIGNED_DATA:
      return "CUFFT_UNALIGNED_DATA";
  }
  return "UNKNOWN";
David Rohr's avatar
David Rohr committed
94 95
}

Luka Stanisic's avatar
Luka Stanisic committed
96 97
/* Handing CUDA Driver errors */

Luka Stanisic's avatar
Luka Stanisic committed
98 99 100 101 102 103 104 105 106 107
#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__;                                                        \
    }                                                                          \
Luka Stanisic's avatar
Luka Stanisic committed
108 109
  } while (false)

Luka Stanisic's avatar
Luka Stanisic committed
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
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";
Luka Stanisic's avatar
Luka Stanisic committed
210 211 212
  }
}

213 214
bioem_cuda::bioem_cuda()
{
Luka Stanisic's avatar
Luka Stanisic committed
215 216 217 218 219 220 221 222
  deviceInitialized = 0;
  GPUAsync = getenv("GPUASYNC") == NULL ? 1 : atoi(getenv("GPUASYNC"));
  GPUWorkload =
      getenv("GPUWORKLOAD") == NULL ? 100 : atoi(getenv("GPUWORKLOAD"));
  if (GPUWorkload == -1)
    GPUWorkload = 100;
  GPUDualStream =
      getenv("GPUDUALSTREAM") == NULL ? 1 : atoi(getenv("GPUDUALSTREAM"));
223 224
}

Luka Stanisic's avatar
Luka Stanisic committed
225
bioem_cuda::~bioem_cuda() { deviceExit(); }
226

Luka Stanisic's avatar
Luka Stanisic committed
227 228 229 230
__global__ void multComplexMap(const mycomplex_t *convmap,
                               const mycomplex_t *refmap, mycuComplex_t *out,
                               const int MapSize, const int maxParallelConv,
                               const int NumberRefMaps, const int Offset)
231
{
Luka Stanisic's avatar
Luka Stanisic committed
232 233 234 235 236 237 238 239 240 241 242 243 244 245 246
  int myConv = myBlockIdxX / NumberRefMaps;
  int myRef = myBlockIdxX - myConv * NumberRefMaps + Offset;
  const mycuComplex_t *myin = (mycuComplex_t *) &refmap[myRef * MapSize];
  const mycuComplex_t *myconv = (mycuComplex_t *) &convmap[myConv * MapSize];
  mycuComplex_t *myout = &out[myBlockIdxX * MapSize];
  for (int i = myThreadIdxX; i < MapSize; i += myBlockDimX)
  {
    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;
  }
247 248
}

Luka Stanisic's avatar
Luka Stanisic committed
249 250 251 252 253
__global__ void
cuDoRefMapsFFT(const int iOrient, const int iConv, const myfloat_t *lCC,
               const myparam5_t *comp_params, bioem_Probability pProb,
               const bioem_param_device param, const bioem_RefMap RefMap,
               const int maxRef, const int Offset)
254
{
Luka Stanisic's avatar
Luka Stanisic committed
255 256 257 258 259 260 261 262
  if (myBlockIdxX * myBlockDimX + myThreadIdxX >= maxRef)
    return;
  const int iRefMap = myBlockIdxX * myBlockDimX + myThreadIdxX + Offset;
  const myfloat_t *mylCC = &lCC[(myBlockIdxX * myBlockDimX + myThreadIdxX) *
                                param.NumberPixels * param.NumberPixels];
  doRefMapFFT(iRefMap, iOrient, iConv, comp_params->amp, comp_params->pha,
              comp_params->env, comp_params->sumC, comp_params->sumsquareC,
              mylCC, pProb, param, RefMap);
263 264
}

Luka Stanisic's avatar
Luka Stanisic committed
265 266 267 268 269 270 271
__global__ void
doRefMap_GPU_Parallel(const int iRefMap, const int iOrient, const int iConv,
                      const int maxParallelConv, const myfloat_t *lCC,
                      const myparam5_t *comp_params, myblockGPU_t *comp_block,
                      bioem_Probability pProb, const bioem_param_device param,
                      const bioem_RefMap RefMap, const int maxRef,
                      const int dispC)
272
{
Luka Stanisic's avatar
Luka Stanisic committed
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
  int myGlobalId = myBlockIdxX * myBlockDimX + myThreadIdxX;
  if (myGlobalId >= maxParallelConv * param.NtotDisp)
    return;
  int myConv = myGlobalId / param.NtotDisp;
  myGlobalId -= myConv * param.NtotDisp;
  int myX = myGlobalId / param.NxDisp;
  myGlobalId -= myX * param.NxDisp;
  int myY = myGlobalId;
  myGlobalId = myBlockIdxX * myBlockDimX + myThreadIdxX;

  int cent_x = (myX * param.GridSpaceCenter + dispC) % param.NumberPixels;
  int cent_y = (myY * param.GridSpaceCenter + dispC) % param.NumberPixels;
  int address = (myConv * maxRef * param.NumberPixels * param.NumberPixels) +
                (cent_x * param.NumberPixels + cent_y);
  myfloat_t value = (myfloat_t) lCC[address] /
                    (myfloat_t)(param.NumberPixels * param.NumberPixels);

  __shared__ myprob_t bestLogpro[CUDA_THREAD_MAX];
  __shared__ int bestId[CUDA_THREAD_MAX];
  __shared__ myprob_t sumExp[CUDA_THREAD_MAX];
  __shared__ myprob_t sumAngles[CUDA_THREAD_MAX];

  int nTotalThreads =
      ((maxParallelConv * param.NtotDisp) < ((myBlockIdxX + 1) * myBlockDimX)) ?
          ((maxParallelConv * param.NtotDisp) - (myBlockIdxX * myBlockDimX)) :
          myBlockDimX;
  int halfPoint = (nTotalThreads + 1) >> 1; // divide by two

  bioem_Probability_map &pProbMap = pProb.getProbMap(iRefMap);

  bestLogpro[myThreadIdxX] =
      calc_logpro(param, comp_params[myConv].amp, comp_params[myConv].pha,
                  comp_params[myConv].env, comp_params[myConv].sumC,
                  comp_params[myConv].sumsquareC, value,
                  RefMap.sum_RefMap[iRefMap], RefMap.sumsquare_RefMap[iRefMap]);
#ifdef DEBUG_PROB
  printf("\t\t\tProb: iRefMap %d, iOrient %d, iConv %d, "
         "cent_x %d, cent_y %d, address %d, value %f, logpro %f\n",
         iRefMap, iOrient, iConv, cent_x, cent_y, address, value,
         bestLogpro[myThreadIdxX]);
#endif
  bestId[myThreadIdxX] = myGlobalId;
  sumExp[myThreadIdxX] = exp(bestLogpro[myThreadIdxX] - pProbMap.Constoadd);
  if (param.writeAngles)
  {
    bioem_Probability_angle &pProbAngle = pProb.getProbAngle(iRefMap, iOrient);
    sumAngles[myThreadIdxX] =
        exp(bestLogpro[myThreadIdxX] - pProbAngle.ConstAngle);
  }
  __syncthreads();
323

Luka Stanisic's avatar
Luka Stanisic committed
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
  // Total number of active threads
  while (nTotalThreads > 1)
  {
    if (myThreadIdxX < (nTotalThreads >> 1))
    {
      // Get the shared value stored by another thread
      myprob_t temp = bestLogpro[myThreadIdxX + halfPoint];
      if (temp > bestLogpro[myThreadIdxX])
      {
        bestLogpro[myThreadIdxX] = temp;
        bestId[myThreadIdxX] = bestId[myThreadIdxX + halfPoint];
      }
      sumExp[myThreadIdxX] += sumExp[myThreadIdxX + halfPoint];
      if (param.writeAngles)
      {
        sumAngles[myThreadIdxX] += sumAngles[myThreadIdxX + halfPoint];
      }
    }
    __syncthreads();
    nTotalThreads = halfPoint;            // divide by two.
    halfPoint = (nTotalThreads + 1) >> 1; // divide by two
    // only the first half of the threads will be active.
  }
  if (myThreadIdxX == 0)
  {
    comp_block[myBlockIdxX].logpro = bestLogpro[0];
    comp_block[myBlockIdxX].id = bestId[0];
    comp_block[myBlockIdxX].sumExp = sumExp[0];
    if (param.writeAngles)
    {
      comp_block[myBlockIdxX].sumAngles = sumAngles[0];
    }
#ifdef DEBUG_PROB
    printf("\t\t\tProb block: iRefMap %d, iOrient %d, iConv %d, "
           "bestlogpro %f, bestId %d, sumExp %f\n",
           iRefMap, iOrient, iConv, bestLogpro[0], bestId[0], sumExp[0]);
#endif
  }
362 363
}

Luka Stanisic's avatar
Luka Stanisic committed
364 365 366 367 368 369 370
__global__ void
doRefMap_GPU_Reduce(const int iRefMap, const int iOrient, const int iConv,
                    const int maxParallelConv, const myfloat_t *lCC,
                    const myparam5_t *comp_params,
                    const myblockGPU_t *comp_block, bioem_Probability pProb,
                    const bioem_param_device param, const bioem_RefMap RefMap,
                    const int maxRef, const int dispC)
371
{
Luka Stanisic's avatar
Luka Stanisic committed
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 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485

  __shared__ myprob_t bestLogpro[CUDA_THREAD_MAX];
  __shared__ int bestId[CUDA_THREAD_MAX];
  __shared__ myprob_t sumExp[CUDA_THREAD_MAX];
  __shared__ myprob_t sumAngles[CUDA_THREAD_MAX];

  // if it is the last block
  int nTotalThreads = myBlockDimX;
  int halfPoint = (nTotalThreads + 1) >> 1; // divide by two

  bioem_Probability_map &pProbMap = pProb.getProbMap(iRefMap);

  bestLogpro[myThreadIdxX] = comp_block[myThreadIdxX].logpro;
  bestId[myThreadIdxX] = comp_block[myThreadIdxX].id;
  sumExp[myThreadIdxX] = comp_block[myThreadIdxX].sumExp;
  if (param.writeAngles)
  {
    sumAngles[myThreadIdxX] = comp_block[myThreadIdxX].sumAngles;
  }
  __syncthreads();
  while (nTotalThreads > 1)
  {
    if (myThreadIdxX < (nTotalThreads >> 1))
    {
      // Get the shared value stored by another thread
      myfloat_t temp = bestLogpro[myThreadIdxX + halfPoint];
      if (temp > bestLogpro[myThreadIdxX])
      {
        bestLogpro[myThreadIdxX] = temp;
        bestId[myThreadIdxX] = bestId[myThreadIdxX + halfPoint];
      }
      sumExp[myThreadIdxX] += sumExp[myThreadIdxX + halfPoint];
      if (param.writeAngles)
      {
        sumAngles[myThreadIdxX] += sumAngles[myThreadIdxX + halfPoint];
      }
    }
    __syncthreads();
    nTotalThreads = halfPoint;            // divide by two.
    halfPoint = (nTotalThreads + 1) >> 1; // divide by two
    // only the first half of the threads will be active.
  }

  if (myThreadIdxX == 0)
  {
    pProbMap.Total += sumExp[0];
    if (pProbMap.Constoadd < bestLogpro[0])
    {
      pProbMap.Total *= exp(-bestLogpro[0] + pProbMap.Constoadd);
      pProbMap.Constoadd = bestLogpro[0];

      // ********** Getting parameters that maximize the probability ***********
      int myGlobalId = bestId[0];
      int myConv = myGlobalId / param.NtotDisp;
      myGlobalId -= myConv * param.NtotDisp;
      int myX = myGlobalId / param.NxDisp;
      myGlobalId -= myX * param.NxDisp;
      int myY = myGlobalId;

      int cent_x = (myX * param.GridSpaceCenter + dispC) % param.NumberPixels;
      int cent_y = (myY * param.GridSpaceCenter + dispC) % param.NumberPixels;
      int address =
          (myConv * maxRef * param.NumberPixels * param.NumberPixels) +
          (cent_x * param.NumberPixels + cent_y);
      myfloat_t value = (myfloat_t) lCC[address] /
                        (myfloat_t)(param.NumberPixels * param.NumberPixels);

      pProbMap.max.max_prob_cent_x =
          -((myX * param.GridSpaceCenter + dispC) - param.NumberPixels);
      pProbMap.max.max_prob_cent_y =
          -((myY * param.GridSpaceCenter + dispC) - param.NumberPixels);
      pProbMap.max.max_prob_orient = iOrient;
      pProbMap.max.max_prob_conv = iConv + myConv;
      pProbMap.max.max_prob_norm =
          -(-comp_params[myConv].sumC * RefMap.sum_RefMap[iRefMap] +
            param.Ntotpi * value) /
          (comp_params[myConv].sumC * comp_params[myConv].sumC -
           comp_params[myConv].sumsquareC * param.Ntotpi);
      pProbMap.max.max_prob_mu =
          -(-comp_params[myConv].sumC * value +
            comp_params[myConv].sumsquareC * RefMap.sum_RefMap[iRefMap]) /
          (comp_params[myConv].sumC * comp_params[myConv].sumC -
           comp_params[myConv].sumsquareC * param.Ntotpi);

#ifdef DEBUG_PROB
      printf("\tProbabilities change: iRefMap %d, iOrient %d, iConv %d, "
             "Total %f, Const %f, bestlogpro %f, sumExp %f, bestId %d\n",
             iRefMap, iOrient, iConv + myConv, pProbMap.Total,
             pProbMap.Constoadd, bestLogpro[0], sumExp[0], bestId[0]);
      printf("\tParameters: myConv %d, myX %d, myY %d, cent_x %d, cent_y %d, "
             "probX %d, probY %d\n",
             myConv, myX, myY, cent_x, cent_y, pProbMap.max.max_prob_cent_x,
             pProbMap.max.max_prob_cent_y);
#endif
    }
#ifdef DEBUG_PROB
    printf("\t\tProbabilities after Reduce: iRefMap %d, iOrient %d, iConv "
           "%d, Total %f, Const %f, bestlogpro %f, sumExp %f, bestId %d\n",
           iRefMap, iOrient, iConv, pProbMap.Total, pProbMap.Constoadd,
           bestLogpro[0], sumExp[0], bestId[0]);
#endif

    if (param.writeAngles)
    {
      bioem_Probability_angle &pProbAngle =
          pProb.getProbAngle(iRefMap, iOrient);
      pProbAngle.forAngles += sumAngles[0];
      if (pProbAngle.ConstAngle < bestLogpro[0])
      {
        pProbAngle.forAngles *= exp(-bestLogpro[0] + pProbAngle.ConstAngle);
        pProbAngle.ConstAngle = bestLogpro[0];
      }
    }
  }
486 487
}

Luka Stanisic's avatar
Luka Stanisic committed
488 489 490 491 492
__global__ void
init_Constoadd(const int iRefMap, const int iOrient, const myfloat_t *lCC,
               const myparam5_t *comp_params, bioem_Probability pProb,
               const bioem_param_device param, const bioem_RefMap RefMap,
               const int initialized_const)
493
{
Luka Stanisic's avatar
Luka Stanisic committed
494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519
  myfloat_t value =
      (myfloat_t) lCC[0] / (myfloat_t)(param.NumberPixels * param.NumberPixels);

  myfloat_t logpro =
      calc_logpro(param, comp_params->amp, comp_params->pha, comp_params->env,
                  comp_params->sumC, comp_params->sumsquareC, value,
                  RefMap.sum_RefMap[iRefMap], RefMap.sumsquare_RefMap[iRefMap]);

  bioem_Probability_map &pProbMap = pProb.getProbMap(iRefMap);

  // Needed only once, in the first projection
  if (!initialized_const)
  {
    pProbMap.Constoadd = logpro;
  }
  // Needed for every projection
  if (param.writeAngles)
  {
    bioem_Probability_angle &pProbAngle = pProb.getProbAngle(iRefMap, iOrient);
    pProbAngle.ConstAngle = logpro;
  }

#ifdef DEBUG_GPU
  printf("\tInitialized pProbMap.Constoadd of refmap %d to %f\n", iRefMap,
         pProbMap.Constoadd);
#endif
520 521
}

Luka Stanisic's avatar
Luka Stanisic committed
522
template <class T> static inline T divup(T num, T divider)
523
{
Luka Stanisic's avatar
Luka Stanisic committed
524
  return ((num + divider - 1) / divider);
525 526
}

Luka Stanisic's avatar
Luka Stanisic committed
527 528 529
int bioem_cuda::compareRefMaps(int iPipeline, int iOrient, int iConv,
                               int maxParallelConv, mycomplex_t *conv_mapsFFT,
                               myparam5_t *comp_params, const int startMap)
530
{
Luka Stanisic's avatar
Luka Stanisic committed
531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682
  if (startMap)
  {
    cout << "Error startMap not implemented for GPU Code\n";
    exit(1);
  }
  printCudaDebugStart();
  if (GPUAsync)
  {
    checkCudaErrors(cudaEventSynchronize(cudaEvent[iPipeline & 1]));
    printCudaDebug("time to synch projections");
  }

  int k = (iPipeline & 1) * param.nTotParallelConv;
  memcpy(&pConvMapFFT_Host[k * param.FFTMapSize],
         conv_mapsFFT[k * param.FFTMapSize],
         param.FFTMapSize * maxParallelConv * sizeof(mycomplex_t));
  printCudaDebug("time for memcpy");
  checkCudaErrors(
      cudaMemcpyAsync(&pConvMapFFT[k * param.FFTMapSize],
                      &pConvMapFFT_Host[k * param.FFTMapSize],
                      param.FFTMapSize * maxParallelConv * sizeof(mycomplex_t),
                      cudaMemcpyHostToDevice, cudaStream[GPUAsync ? 2 : 0]));
  // If one wants just a single tranfer, without memcpy:
  // checkCudaErrors(cudaMemcpyAsync(&pConvMapFFT[k * param.FFTMapSize],
  // conv_mapsFFT[k * param.FFTMapSize], param.FFTMapSize * maxParallelConv *
  // sizeof(mycomplex_t), cudaMemcpyHostToDevice, cudaStream[GPUAsync ? 2 :
  // 0]));
  checkCudaErrors(cudaMemcpyAsync(&pTmp_comp_params[k], &comp_params[k],
                                  maxParallelConv * sizeof(myparam5_t),
                                  cudaMemcpyHostToDevice,
                                  cudaStream[GPUAsync ? 2 : 0]));
  printCudaDebug("time for asyncmemcpy");
  if (GPUAsync)
  {
    checkCudaErrors(cudaEventRecord(cudaEvent[2], cudaStream[2]));
    checkCudaErrors(cudaStreamWaitEvent(cudaStream[0], cudaEvent[2], 0));
  }
  if (GPUDualStream)
  {
    checkCudaErrors(cudaEventRecord(cudaFFTEvent[0], cudaStream[0]));
    checkCudaErrors(cudaStreamWaitEvent(cudaStream[1], cudaFFTEvent[0], 0));
  }
  for (int offset = 0, stream = 0; offset < maxRef;
       offset += param.nTotParallelMaps, stream++)
  {
    if (!GPUDualStream)
      stream = 0;
    const int nRef = min(param.nTotParallelMaps, maxRef - offset);
    multComplexMap<<<maxParallelConv * nRef, CudaThreadCount, 0,
                     cudaStream[stream & 1]>>>(
        &pConvMapFFT[k * param.FFTMapSize], pRefMapsFFT, pFFTtmp2[stream & 1],
        param.FFTMapSize, maxParallelConv, nRef, offset);
    printCudaDebug("time for multComplexMap kernel");
    cufftResult err = mycufftExecC2R(offset + param.nTotParallelMaps > maxRef ?
                                         plan[1][stream & 1] :
                                         plan[0][stream & 1],
                                     pFFTtmp2[stream & 1], pFFTtmp[stream & 1]);
    if (err != CUFFT_SUCCESS)
    {
      cout << "Error running CUFFT " << cufftGetErrorStrung(err) << "\n";
      exit(1);
    }
    printCudaDebug("time for mycufftExecC2R kernel");
    if (BioEMAlgo == 1)
    {
      for (int conv = 0; conv < maxParallelConv; conv++)
      {
        cuDoRefMapsFFT<<<divup(nRef, CudaThreadCount), CudaThreadCount, 0,
                         cudaStream[stream & 1]>>>(
            iOrient, iConv + conv,
            pFFTtmp[stream & 1] +
                conv * nRef * param.param_device.NumberPixels *
                    param.param_device.NumberPixels,
            &pTmp_comp_params[k + conv], pProb_device, param.param_device,
            *gpumap, nRef, offset);
        printCudaDebug("time for cuDoRefMapsFFT kernel");
      }
    }
    else
    {
      for (int refmap = offset; refmap < nRef + offset; refmap++)
      {
        // First iteration needs to initialize Constoadd with the first valid
        // value to avoid overflow due to high sumExp values
        if ((initialized_const[refmap] == false) ||
            (param.param_device.writeAngles && iConv == 0))
        {
          init_Constoadd<<<1, 1, 0, cudaStream[stream & 1]>>>(
              refmap, iOrient,
              pFFTtmp[stream & 1] +
                  (refmap - offset) * param.param_device.NumberPixels *
                      param.param_device.NumberPixels,
              &pTmp_comp_params[k], pProb_device, param.param_device, *gpumap,
              (int) initialized_const[refmap]);
          initialized_const[refmap] = true;
          printCudaDebug("time for init_Constoadd kernel");
        }

        doRefMap_GPU_Parallel<<<divup(maxParallelConv *
                                          param.param_device.NtotDisp,
                                      CudaThreadCount),
                                CudaThreadCount, 0, cudaStream[stream & 1]>>>(
            refmap, iOrient, iConv, maxParallelConv,
            pFFTtmp[stream & 1] +
                (refmap - offset) * param.param_device.NumberPixels *
                    param.param_device.NumberPixels,
            &pTmp_comp_params[k], &pTmp_comp_blocks[refmap * Ncomp_blocks],
            pProb_device, param.param_device, *gpumap, nRef,
            param.param_device.NumberPixels -
                param.param_device.maxDisplaceCenter);
        printCudaDebug("time for doRefMaps_GPU_Parallel kernel");

        doRefMap_GPU_Reduce<<<1, divup(maxParallelConv *
                                           param.param_device.NtotDisp,
                                       CudaThreadCount),
                              0, cudaStream[stream & 1]>>>(
            refmap, iOrient, iConv, maxParallelConv,
            pFFTtmp[stream & 1] +
                (refmap - offset) * param.param_device.NumberPixels *
                    param.param_device.NumberPixels,
            &pTmp_comp_params[k], &pTmp_comp_blocks[refmap * Ncomp_blocks],
            pProb_device, param.param_device, *gpumap, nRef,
            param.param_device.NumberPixels -
                param.param_device.maxDisplaceCenter);
        printCudaDebug("time for doRefMaps_GPU_Reduce kernel");
      }
    }
  }
  checkCudaErrors(cudaPeekAtLastError());

  if (GPUDualStream)
  {
    checkCudaErrors(cudaEventRecord(cudaFFTEvent[1], cudaStream[1]));
    checkCudaErrors(cudaStreamWaitEvent(cudaStream[0], cudaFFTEvent[1], 0));
  }

  if ((BioEMAlgo == 1) && (GPUWorkload < 100))
  {
    bioem::compareRefMaps(iPipeline, iOrient, iConv, maxParallelConv,
                          conv_mapsFFT, comp_params, maxRef);
    printCudaDebug("time to run OMP");
  }
  if (GPUAsync)
  {
    checkCudaErrors(cudaEventRecord(cudaEvent[iPipeline & 1], cudaStream[0]));
  }
  else
  {
    checkCudaErrors(cudaStreamSynchronize(cudaStream[0]));
    printCudaDebug("time to synch at the end");
  }
  return (0);
683 684
}

David Rohr's avatar
David Rohr committed
685 686
int bioem_cuda::selectCudaDevice()
{
Luka Stanisic's avatar
Luka Stanisic committed
687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725
  int count;
  int bestDevice = 0;
  cudaDeviceProp deviceProp;

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

  /* Get number of available CUDA devices */
  checkCudaErrors(cudaGetDeviceCount(&count));
  if (count == 0)
  {
    printf("No CUDA device detected\n");
    return (1);
  }

  /* Find the best GPU */
  long long int bestDeviceSpeed = -1, deviceSpeed = -1;
  for (int i = 0; i < count; i++)
  {
    cudaGetDeviceProperties(&deviceProp, i);
    deviceSpeed = (long long int) deviceProp.multiProcessorCount *
                  (long long int) deviceProp.clockRate *
                  (long long int) deviceProp.warpSize;
    if (deviceSpeed > bestDeviceSpeed)
    {
      bestDevice = i;
      bestDeviceSpeed = deviceSpeed;
    }
  }

  /* Get user-specified GPU choice */
  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);
    }
726
#ifdef WITH_MPI
Luka Stanisic's avatar
Luka Stanisic committed
727 728 729 730
    if (device == -1)
    {
      device = mpi_rank % count;
    }
731
#endif
Luka Stanisic's avatar
Luka Stanisic committed
732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786
    if (device < 0)
    {
      printf("Negative CUDA device specified: %d, invalid!\n", device);
      exit(1);
    }
    bestDevice = device;
  }

  /* 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());
  }

  /* Debugging information about CUDA devices used by the current process */
  if (DebugOutput >= 2)
  {
    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);
    printf("computeMode = %d\n", deviceProp.computeMode);
787
#if CUDA_VERSION > 3010
Luka Stanisic's avatar
Luka Stanisic committed
788
    size_t free, total;
789
#else
Luka Stanisic's avatar
Luka Stanisic committed
790
    unsigned int free, total;
791
#endif
Luka Stanisic's avatar
Luka Stanisic committed
792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815
    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);
  }

  if (DebugOutput >= 1)
  {
    printf("BioEM for CUDA initialized (MPI Rank %d), %d GPUs found, using GPU "
           "%d\n",
           mpi_rank, count, bestDevice);
  }

  return (0);
David Rohr's avatar
David Rohr committed
816 817
}

818 819
int bioem_cuda::deviceInit()
{
Luka Stanisic's avatar
Luka Stanisic committed
820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912
  deviceExit();

  selectCudaDevice();

  gpumap = new bioem_RefMap;
  memcpy(gpumap, &RefMap, sizeof(bioem_RefMap));

  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->sum_RefMap = sum;
  gpumap->sumsquare_RefMap = sumsquare;

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

  for (int i = 0; i < PIPELINE_LVL; i++)
  {
    checkCudaErrors(cudaStreamCreate(&cudaStream[i]));
    checkCudaErrors(cudaEventCreate(&cudaEvent[i]));
  }
  for (int i = 0; i < MULTISTREAM_LVL; i++)
  {
    checkCudaErrors(cudaEventCreate(&cudaFFTEvent[i]));
  }
  if (GPUAsync)
  {
    checkCudaErrors(cudaStreamCreate(&cudaStream[2]));
    checkCudaErrors(cudaEventCreate(&cudaEvent[2]));
  }

  checkCudaErrors(
      cudaMalloc(&pRefMapsFFT,
                 RefMap.ntotRefMap * param.FFTMapSize * sizeof(mycomplex_t)));
  checkCudaErrors(
      cudaMalloc(&pFFTtmp2[0], param.nTotParallelConv * param.nTotParallelMaps *
                                   param.FFTMapSize * MULTISTREAM_LVL *
                                   sizeof(mycomplex_t)));
  checkCudaErrors(
      cudaMalloc(&pFFTtmp[0], param.nTotParallelConv * param.nTotParallelMaps *
                                  param.param_device.NumberPixels *
                                  param.param_device.NumberPixels *
                                  MULTISTREAM_LVL * sizeof(myfloat_t)));
  for (int i = 1; i < MULTISTREAM_LVL; i++)
  {
    pFFTtmp2[i] =
        pFFTtmp2[0] +
        i * param.nTotParallelConv * param.nTotParallelMaps * param.FFTMapSize;
    pFFTtmp[i] = pFFTtmp[0] +
                 i * param.nTotParallelConv * param.nTotParallelMaps *
                     param.param_device.NumberPixels *
                     param.param_device.NumberPixels;
  }
  checkCudaErrors(cudaMalloc(&pConvMapFFT, param.nTotParallelConv *
                                               param.FFTMapSize * PIPELINE_LVL *
                                               sizeof(mycomplex_t)));
  checkCudaErrors(cudaHostAlloc(&pConvMapFFT_Host,
                                param.nTotParallelConv * param.FFTMapSize *
                                    PIPELINE_LVL * sizeof(mycomplex_t),
                                0));
  checkCudaErrors(
      cudaMemcpy(pRefMapsFFT, RefMap.RefMapsFFT,
                 RefMap.ntotRefMap * param.FFTMapSize * sizeof(mycomplex_t),
                 cudaMemcpyHostToDevice));
  checkCudaErrors(
      cudaMalloc(&pTmp_comp_params,
                 param.nTotParallelConv * PIPELINE_LVL * sizeof(myparam5_t)));
  Ncomp_blocks = divup(param.nTotParallelConv * param.param_device.NtotDisp,
                       CudaThreadCount);
  if (Ncomp_blocks > CudaThreadCount)
  {
    cout << "Error with input parameters. Check CudaThreadCount, "
            "displacements and max number of parallel comparisons\n";
    exit(1);
  }
  checkCudaErrors(
      cudaMalloc(&pTmp_comp_blocks,
                 Ncomp_blocks * RefMap.ntotRefMap * sizeof(myblockGPU_t)));

  initialized_const = new bool[RefMap.ntotRefMap];
  for (int i = 0; i < RefMap.ntotRefMap; i++)
    initialized_const[i] = false;

  deviceInitialized = 1;
  return (0);
913 914 915 916
}

int bioem_cuda::deviceExit()
{
Luka Stanisic's avatar
Luka Stanisic committed
917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952
  if (deviceInitialized == 0)
    return (0);

  cudaFree(pProb_memory);
  cudaFree(sum);
  cudaFree(sumsquare);
  for (int i = 0; i < PIPELINE_LVL; i++)
  {
    cudaStreamDestroy(cudaStream[i]);
    cudaEventDestroy(cudaEvent[i]);
  }
  for (int i = 0; i < MULTISTREAM_LVL; i++)
  {
    cudaEventDestroy(cudaFFTEvent[i]);
  }

  cudaFree(pRefMapsFFT);
  cudaFree(pConvMapFFT);
  cudaFreeHost(pConvMapFFT_Host);
  cudaFree(pFFTtmp[0]);
  cudaFree(pFFTtmp2[0]);
  cudaFree(pTmp_comp_params);
  cudaFree(pTmp_comp_blocks);

  if (GPUAsync)
  {
    cudaStreamDestroy(cudaStream[2]);
    cudaEventDestroy(cudaEvent[2]);
  }

  delete gpumap;
  delete initialized_const;
  cudaThreadExit();

  deviceInitialized = 0;
  return (0);
953 954 955 956
}

int bioem_cuda::deviceStartRun()
{
Luka Stanisic's avatar
Luka Stanisic committed
957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021
  if (GPUWorkload >= 100)
  {
    maxRef = RefMap.ntotRefMap;
    pProb_host = &pProb;
  }
  else
  {
    maxRef = ((size_t) RefMap.ntotRefMap * (size_t) GPUWorkload / 100) < 1 ?
                 (size_t) RefMap.ntotRefMap :
                 (size_t) RefMap.ntotRefMap * (size_t) GPUWorkload / 100;
    pProb_host = new bioem_Probability;
    pProb_host->init(maxRef, param.nTotGridAngles, *this);
    pProb_host->copyFrom(&pProb, *this);
  }

  pProb_device = *pProb_host;
  pProb_device.ptr = pProb_memory;
  pProb_device.set_pointers();
  checkCudaErrors(
      cudaMemcpyAsync(pProb_device.ptr, pProb_host->ptr,
                      pProb_host->get_size(maxRef, param.nTotGridAngles,
                                           param.param_device.writeAngles),
                      cudaMemcpyHostToDevice, cudaStream[0]));

  if (maxRef / (param.nTotParallelMaps * param.nTotParallelConv) >
      (double) SPLIT_MAPS_LVL)
  {
    cout << "Error planning CUFFT dimensions\n";
    exit(1);
  }
  for (int j = 0; j < MULTISTREAM_LVL; j++)
  {
    for (int i = 0; i < SPLIT_MAPS_LVL; i++)
    {
      if (i && maxRef % param.nTotParallelMaps == 0)
        continue;
      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 % param.nTotParallelMaps) * param.nTotParallelConv) :
                  (param.nTotParallelMaps * param.nTotParallelConv)) !=
          CUFFT_SUCCESS)
      {
        cout << "Error planning CUFFT\n";
        exit(1);
      }
      if (cufftSetCompatibilityMode(
              plan[i][j], CUFFT_COMPATIBILITY_FFTW_PADDING) != 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);
      }
    }
    if (!GPUDualStream)
      break;
  }

  return (0);
1022 1023 1024 1025
}

int bioem_cuda::deviceFinishRun()
{
Luka Stanisic's avatar
Luka Stanisic committed
1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044
  if (GPUAsync)
    cudaStreamSynchronize(cudaStream[0]);
  checkCudaErrors(
      cudaMemcpyAsync(pProb_host->ptr, pProb_device.ptr,
                      pProb_host->get_size(maxRef, param.nTotGridAngles,
                                           param.param_device.writeAngles),
                      cudaMemcpyDeviceToHost, cudaStream[0]));

  for (int j = 0; j < MULTISTREAM_LVL; j++)
  {
    for (int i = 0; i < SPLIT_MAPS_LVL; i++)
    {
      if (i && maxRef % param.nTotParallelMaps == 0)
        continue;
      cufftDestroy(plan[i][j]);
    }
    if (!GPUDualStream)
      break;
  }
1045

Luka Stanisic's avatar
Luka Stanisic committed
1046 1047 1048 1049 1050 1051 1052 1053 1054
  cudaThreadSynchronize();
  if (GPUWorkload < 100)
  {
    pProb.copyFrom(pProb_host, *this);
    free_device_host(pProb_host->ptr);
    delete[] pProb_host;
  }

  return (0);
1055 1056
}

Luka Stanisic's avatar
Luka Stanisic committed
1057
void *bioem_cuda::malloc_device_host(size_t size)
1058
{
Luka Stanisic's avatar
Luka Stanisic committed
1059 1060 1061
  void *ptr;
  checkCudaErrors(cudaHostAlloc(&ptr, size, 0));
  return (ptr);
1062 1063
}

Luka Stanisic's avatar
Luka Stanisic committed
1064 1065
void bioem_cuda::free_device_host(void *ptr) { cudaFreeHost(ptr); }

Luka Stanisic's avatar
Luka Stanisic committed
1066 1067
void bioem_cuda::rebalance(int workload)
{
Luka Stanisic's avatar
Luka Stanisic committed
1068 1069
  if ((workload < 0) || (workload > 100) || (workload == GPUWorkload))
    return;
Luka Stanisic's avatar
Luka Stanisic committed
1070

Luka Stanisic's avatar
Luka Stanisic committed
1071
  deviceFinishRun();
Luka Stanisic's avatar
Luka Stanisic committed
1072

Luka Stanisic's avatar
Luka Stanisic committed
1073 1074 1075 1076
  if (DebugOutput >= 2)
  {
    printf("\t\tSetting GPU workload to %d%% (rank %d)\n", workload, mpi_rank);
  }
Luka Stanisic's avatar
Luka Stanisic committed
1077

Luka Stanisic's avatar
Luka Stanisic committed
1078 1079
  GPUWorkload = workload;
  maxRef = (size_t) RefMap.ntotRefMap * (size_t) GPUWorkload / 100;
Luka Stanisic's avatar
Luka Stanisic committed
1080

Luka Stanisic's avatar
Luka Stanisic committed
1081
  deviceStartRun();
Luka Stanisic's avatar
Luka Stanisic committed
1082 1083
}

Luka Stanisic's avatar
Luka Stanisic committed
1084
bioem *bioem_cuda_create()
1085
{
Luka Stanisic's avatar
Luka Stanisic committed
1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096
  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;
  }

  return new bioem_cuda;
1097
}