InstanceCUDA.cpp 44.8 KB
Newer Older
1 2 3
// Copyright (C) 2020 ASTRON (Netherlands Institute for Radio Astronomy)
// SPDX-License-Identifier: GPL-3.0-or-later

4
#include <iomanip>  // setprecision
5

6
#include "InstanceCUDA.h"
Bram Veenboer's avatar
Bram Veenboer committed
7
#include "PowerRecord.h"
8

9
using namespace idg::kernel;
10
using namespace powersensor;
11

12 13
#define NR_CORRELATIONS 4

14 15 16 17 18
/*
 * Option to enable repeated kernel invocations
 * this is used to measure energy consumpton
 * using a low-resolution power measurement (NVML)
 */
19 20 21 22
#define ENABLE_REPEAT_KERNELS 0
#define NR_REPETITIONS_GRIDDER 10
#define NR_REPETITIONS_ADDER 50
#define NR_REPETITIONS_GRID_FFT 500
23

24 25 26 27 28
/*
 * Use custom FFT kernel
 */
#define USE_CUSTOM_FFT 0

29
namespace idg {
30 31 32 33 34 35 36 37 38 39 40 41 42
namespace kernel {
namespace cuda {

// Constructor
InstanceCUDA::InstanceCUDA(ProxyInfo& info, int device_nr, int device_id)
    : KernelsInstance(), mInfo(info) {
#if defined(DEBUG)
  std::cout << __func__ << std::endl;
#endif

  // Initialize members
  device.reset(new cu::Device(device_id));
  context.reset(new cu::Context(*device));
43 44 45
  executestream.reset(new cu::Stream(*context));
  htodstream.reset(new cu::Stream(*context));
  dtohstream.reset(new cu::Stream(*context));
46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61

  // Set kernel parameters
  set_parameters();

  // Compile kernels
  compile_kernels();

  // Load kernels
  load_kernels();

  // Initialize power sensor
  powerSensor = get_power_sensor(sensor_device, device_nr);
}

// Destructor
InstanceCUDA::~InstanceCUDA() {
62 63 64 65
#if defined(DEBUG)
  std::cout << __func__ << std::endl;
#endif

66 67 68
  free_host_memory();
  free_device_memory();
  free_fft_plans();
69
  free_events();
70 71 72 73 74
  mModules.clear();
  executestream.reset();
  htodstream.reset();
  dtohstream.reset();
  context.reset();
75
  device.reset();
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
  delete powerSensor;
}

/*
    Compilation
*/
std::string InstanceCUDA::get_compiler_flags() {
  // Constants
  std::stringstream flags_constants;
  flags_constants << "-DNR_POLARIZATIONS=" << NR_CORRELATIONS;

  // CUDA specific flags
  std::stringstream flags_cuda;
  flags_cuda << "-use_fast_math ";
#if defined(CUDA_KERNEL_DEBUG)
  flags_cuda << " -G ";
#else
  flags_cuda << "-lineinfo ";
#endif
  flags_cuda << "-src-in-ptx";

  // Device specific flags
  int capability = (*device).get_capability();
  std::stringstream flags_device;
  flags_device << "-arch=sm_" << capability;

  // Include flags
  std::stringstream flags_includes;
  flags_includes << "-I" << IDG_INSTALL_DIR << "/include";

  // Combine flags
  std::string flags = " " + flags_cuda.str() + " " + flags_device.str() + " " +
                      flags_constants.str() + " " + flags_includes.str();
  return flags;
}

cu::Module* InstanceCUDA::compile_kernel(std::string& flags, std::string& src,
                                         std::string& bin) {
  // Create a string with the full path to the cubin file "kernel.cubin"
  std::string lib = mInfo.get_path_to_lib() + "/" + bin;

  // Create a string for all sources that are combined
  std::string source = mInfo.get_path_to_src() + "/" + src;

  // Call the compiler
  cu::Source(source.c_str()).compile(lib.c_str(), flags.c_str());

  // Create module
124
  return new cu::Module(*context, lib.c_str());
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
}

void InstanceCUDA::compile_kernels() {
#if defined(DEBUG)
  std::cout << __func__ << std::endl;
#endif

  // Get source directory
  std::string srcdir = auxiliary::get_lib_dir() + "/idg-cuda";
#if defined(DEBUG)
  std::cout << "Searching for source files in: " << srcdir << std::endl;
#endif

  // Create temp directory
  char tmpdir[] = "/tmp/idg-XXXXXX";
  char* tmpdir_ = mkdtemp(tmpdir);
  if (!tmpdir_) {
    throw std::runtime_error("could not create tmp directory.");
  }
#if defined(DEBUG)
  std::cout << "Temporary files will be stored in: " << tmpdir << std::endl;
#endif

  // Get compiler flags
  std::string flags_common = get_compiler_flags();

  // Create vector of source filenames, filenames and flags
  std::vector<std::string> src;
  std::vector<std::string> cubin;
  std::vector<std::string> flags;

  // Gridder
  src.push_back("KernelGridder.cu");
  cubin.push_back("Gridder.cubin");
  std::stringstream flags_gridder;
  flags_gridder << flags_common;
  flags_gridder << " -DBATCH_SIZE=" << batch_gridder;
  flags.push_back(flags_gridder.str());

  // Degridder
  src.push_back("KernelDegridder.cu");
  cubin.push_back("Degridder.cubin");
  std::stringstream flags_degridder;
  flags_degridder << flags_common;
  flags_degridder << " -DBATCH_SIZE=" << batch_degridder;
  flags.push_back(flags_degridder.str());

  // Scaler
  src.push_back("KernelScaler.cu");
  cubin.push_back("Scaler.cubin");
  flags.push_back(flags_common);

  // Adder
  src.push_back("KernelAdder.cu");
  cubin.push_back("Adder.cubin");
  std::stringstream flags_adder;
  flags_adder << flags_common;
  flags_adder << " -DTILE_SIZE_GRID=" << tile_size_grid;
  flags.push_back(flags_adder.str());

  // Splitter
  src.push_back("KernelSplitter.cu");
  cubin.push_back("Splitter.cubin");
  std::stringstream flags_splitter;
  flags_splitter << flags_common;
  flags_splitter << " -DTILE_SIZE_GRID=" << tile_size_grid;
  flags.push_back(flags_splitter.str());

  // Calibrate
  src.push_back("KernelCalibrate.cu");
  cubin.push_back("Calibrate.cubin");
  flags.push_back(flags_common);

198 199 200 201 202
  // Average beam
  src.push_back("KernelAverageBeam.cu");
  cubin.push_back("AverageBeam.cubin");
  flags.push_back(flags_common);

Bram Veenboer's avatar
Bram Veenboer committed
203 204 205 206 207
  // FFT shift
  src.push_back("KernelFFTShift.cu");
  cubin.push_back("KernelFFTShift.cubin");
  flags.push_back(flags_common);

208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231
// FFT
#if USE_CUSTOM_FFT
  src.push_back("KernelFFT.cu");
  cubin.push_back("FFT.cubin");
  flags.push_back(flags_common);
#endif

  // Compile all kernels
  for (unsigned i = 0; i < src.size(); i++) {
    mModules.push_back(std::unique_ptr<cu::Module>());
  }
#pragma omp parallel for
  for (unsigned i = 0; i < src.size(); i++) {
    mModules[i].reset(compile_kernel(flags[i], src[i], cubin[i]));
  }
}

void InstanceCUDA::load_kernels() {
  CUfunction function;
  unsigned found = 0;

  // Load gridder function
  if (cuModuleGetFunction(&function, *mModules[0], name_gridder.c_str()) ==
      CUDA_SUCCESS) {
232
    function_gridder.reset(new cu::Function(*context, function));
233 234 235 236 237 238
    found++;
  }

  // Load degridder function
  if (cuModuleGetFunction(&function, *mModules[1], name_degridder.c_str()) ==
      CUDA_SUCCESS) {
239
    function_degridder.reset(new cu::Function(*context, function));
240 241 242 243 244 245
    found++;
  }

  // Load scalar function
  if (cuModuleGetFunction(&function, *mModules[2], name_scaler.c_str()) ==
      CUDA_SUCCESS) {
246
    function_scaler.reset(new cu::Function(*context, function));
247 248 249 250 251 252
    found++;
  }

  // Load adder function
  if (cuModuleGetFunction(&function, *mModules[3], name_adder.c_str()) ==
      CUDA_SUCCESS) {
253
    function_adder.reset(new cu::Function(*context, function));
254 255 256 257 258 259
    found++;
  }

  // Load splitter function
  if (cuModuleGetFunction(&function, *mModules[4], name_splitter.c_str()) ==
      CUDA_SUCCESS) {
260
    function_splitter.reset(new cu::Function(*context, function));
261 262 263 264 265 266
    found++;
  }

  // Load calibration functions
  if (cuModuleGetFunction(&function, *mModules[5],
                          name_calibrate_lmnp.c_str()) == CUDA_SUCCESS) {
Bram Veenboer's avatar
Bram Veenboer committed
267
    functions_calibrate.emplace_back(new cu::Function(*context, function));
268 269 270 271
    found++;
  }
  if (cuModuleGetFunction(&function, *mModules[5],
                          name_calibrate_sums.c_str()) == CUDA_SUCCESS) {
Bram Veenboer's avatar
Bram Veenboer committed
272
    functions_calibrate.emplace_back(new cu::Function(*context, function));
273 274 275
  }
  if (cuModuleGetFunction(&function, *mModules[5],
                          name_calibrate_gradient.c_str()) == CUDA_SUCCESS) {
Bram Veenboer's avatar
Bram Veenboer committed
276
    functions_calibrate.emplace_back(new cu::Function(*context, function));
277 278 279
  }
  if (cuModuleGetFunction(&function, *mModules[5],
                          name_calibrate_hessian.c_str()) == CUDA_SUCCESS) {
Bram Veenboer's avatar
Bram Veenboer committed
280
    functions_calibrate.emplace_back(new cu::Function(*context, function));
281 282
  }

283 284 285 286 287 288 289
  // Load average beam function
  if (cuModuleGetFunction(&function, *mModules[6], name_average_beam.c_str()) ==
      CUDA_SUCCESS) {
    function_average_beam.reset(new cu::Function(*context, function));
    found++;
  }

Bram Veenboer's avatar
Bram Veenboer committed
290 291 292 293 294 295 296
  // Load FFT shift function
  if (cuModuleGetFunction(&function, *mModules[7], name_fft_shift.c_str()) ==
      CUDA_SUCCESS) {
    function_fft_shift.reset(new cu::Function(*context, function));
    found++;
  }

297 298
// Load FFT function
#if USE_CUSTOM_FFT
Bram Veenboer's avatar
Bram Veenboer committed
299
  if (cuModuleGetFunction(&function, *mModules[8], name_fft.c_str()) ==
300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 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 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388
      CUDA_SUCCESS) {
    function_fft.reset(new cu::Function(function));
    found++;
  }
#endif

  // Verify that all functions are found
  if (found != mModules.size()) {
    std::cerr << "Incorrect number of functions found: " << found
              << " != " << mModules.size() << std::endl;
    exit(EXIT_FAILURE);
  }
}

void InstanceCUDA::set_parameters_kepler() {
  batch_gridder = 192;
  batch_degridder = 192;
}

void InstanceCUDA::set_parameters_maxwell() {
  batch_gridder = 384;
  batch_degridder = 512;
}

void InstanceCUDA::set_parameters_gp100() {
  batch_gridder = 256;
  batch_degridder = 128;
}

void InstanceCUDA::set_parameters_pascal() {
  batch_gridder = 384;
  batch_degridder = 512;
}

void InstanceCUDA::set_parameters_volta() {
  batch_gridder = 128;
  batch_degridder = 256;
}

void InstanceCUDA::set_parameters_default() {
  block_gridder = dim3(128);
  block_degridder = dim3(128);
  block_calibrate = dim3(128);
  block_adder = dim3(128);
  block_splitter = dim3(128);
  block_scaler = dim3(128);
  tile_size_grid = 128;
}

void InstanceCUDA::set_parameters() {
#if defined(DEBUG)
  std::cout << __func__ << std::endl;
#endif

  int capability = (*device).get_capability();

  set_parameters_default();

  if (capability >= 70) {
    set_parameters_volta();
  } else if (capability >= 61) {
    set_parameters_pascal();
  } else if (capability == 60) {
    set_parameters_gp100();
  } else if (capability >= 50) {
    set_parameters_maxwell();
  } else if (capability >= 30) {
    set_parameters_kepler();
  } else {
    // IDG has never been tested with pre-Kepler GPUs
    set_parameters_kepler();
  }

  // Override parameters from environment
  char* cstr_batch_size = getenv("BATCHSIZE");
  if (cstr_batch_size) {
    auto batch_size = atoi(cstr_batch_size);
    batch_gridder = batch_size;
    batch_degridder = batch_size;
  }
  char* cstr_block_size = getenv("BLOCKSIZE");
  if (cstr_block_size) {
    auto block_size = atoi(cstr_block_size);
    block_gridder = dim3(block_size);
    block_degridder = dim3(block_size);
  }
}

std::ostream& operator<<(std::ostream& os, InstanceCUDA& d) {
389
  cu::ScopedContext scc(d.get_context());
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
  os << d.get_device().get_name() << std::endl;
  os << std::setprecision(2);
  os << std::fixed;

  // Device memory
  auto device_memory_total =
      d.get_device().get_total_memory() / (1024 * 1024);  // MBytes
  auto device_memory_free =
      d.get_device().get_free_memory() / (1024 * 1024);  // MBytes
  os << "\tDevice memory : " << device_memory_free << " Mb  / "
     << device_memory_total << " Mb (free / total)" << std::endl;

  // Shared memory
  auto shared_memory =
      d.get_device()
          .get_attribute<
              CU_DEVICE_ATTRIBUTE_MAX_SHARED_MEMORY_PER_BLOCK>();  // Bytes
  os << "\tShared memory : " << shared_memory / (float)1024 << " Kb"
     << std::endl;

  // Frequencies
  auto clock_frequency =
      d.get_device().get_attribute<CU_DEVICE_ATTRIBUTE_CLOCK_RATE>() /
      1000;  // Mhz
  os << "\tClk frequency : " << clock_frequency << " Ghz" << std::endl;
  auto mem_frequency =
      d.get_device().get_attribute<CU_DEVICE_ATTRIBUTE_MEMORY_CLOCK_RATE>() /
      1000;  // Mhz
  os << "\tMem frequency : " << mem_frequency << " Ghz" << std::endl;

  // Cores/bus
  auto nr_sm =
      d.get_device().get_attribute<CU_DEVICE_ATTRIBUTE_MULTIPROCESSOR_COUNT>();
  auto mem_bus_width =
      d.get_device()
          .get_attribute<
              CU_DEVICE_ATTRIBUTE_GLOBAL_MEMORY_BUS_WIDTH>();  // Bits
  os << "\tNumber of SM  : " << nr_sm << std::endl;
  os << "\tMem bus width : " << mem_bus_width << " bit" << std::endl;
  os << "\tMem bandwidth : " << 2 * (mem_bus_width / 8) * mem_frequency / 1000
     << " GB/s" << std::endl;

  auto nr_threads =
      d.get_device()
          .get_attribute<CU_DEVICE_ATTRIBUTE_MAX_THREADS_PER_MULTIPROCESSOR>();
  os << "\tNumber of threads  : " << nr_threads << std::endl;

  // Misc
  os << "\tCapability    : " << d.get_device().get_capability() << std::endl;

  // Unified memory
  auto supports_managed_memory =
      d.get_device().get_attribute<CU_DEVICE_ATTRIBUTE_MANAGED_MEMORY>();
  os << "\tUnified memory : " << supports_managed_memory << std::endl;

  os << std::endl;
  return os;
}

State InstanceCUDA::measure() { return powerSensor->read(); }

void InstanceCUDA::measure(PowerRecord& record, cu::Stream& stream) {
  record.sensor = powerSensor;
  record.enqueue(stream);
}

457
cu::Event& InstanceCUDA::get_event() {
Bram Veenboer's avatar
Bram Veenboer committed
458 459 460 461 462 463
  // Create new event
  cu::Event* event = new cu::Event(*context);

  // This event is used in a callback, where it can not be destroyed
  // after use. Instead, register the event globally, and take care of
  // destruction there.
464
  events.emplace_back(event);
Bram Veenboer's avatar
Bram Veenboer committed
465 466

  // Return a reference to the event
467 468 469
  return *event;
}

470 471 472 473 474 475 476
typedef struct {
  PowerRecord* start;
  PowerRecord* end;
  Report* report;
  void (Report::*update_report)(State&, State&);
} UpdateData;

477 478
UpdateData* get_update_data(cu::Event& event, PowerSensor* sensor,
                            Report* report,
479 480
                            void (Report::*update_report)(State&, State&)) {
  UpdateData* data = new UpdateData();
481 482
  data->start = new PowerRecord(event, sensor);
  data->end = new PowerRecord(event, sensor);
483 484 485 486 487 488 489 490 491 492 493 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 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534
  data->report = report;
  data->update_report = update_report;
  return data;
}

void update_report_callback(CUstream, CUresult, void* userData) {
  UpdateData* data = static_cast<UpdateData*>(userData);
  PowerRecord* start = data->start;
  PowerRecord* end = data->end;
  Report* report = data->report;
  (report->*data->update_report)(start->state, end->state);
  delete start;
  delete end;
  delete data;
}

void InstanceCUDA::start_measurement(void* ptr) {
  UpdateData* data = (UpdateData*)ptr;

  // Schedule the first measurement (prior to kernel execution)
  data->start->enqueue(*executestream);
}

void InstanceCUDA::end_measurement(void* ptr) {
  UpdateData* data = (UpdateData*)ptr;

  // Schedule the second measurement (after the kernel execution)
  data->end->enqueue(*executestream);

  // Afterwards, update the report according to the two measurements
  executestream->addCallback((CUstreamCallback)&update_report_callback, data);
}

void InstanceCUDA::launch_gridder(
    int time_offset, int nr_subgrids, int grid_size, int subgrid_size,
    float image_size, float w_step, int nr_channels, int nr_stations,
    cu::DeviceMemory& d_uvw, cu::DeviceMemory& d_wavenumbers,
    cu::DeviceMemory& d_visibilities, cu::DeviceMemory& d_spheroidal,
    cu::DeviceMemory& d_aterm, cu::DeviceMemory& d_aterm_indices,
    cu::DeviceMemory& d_avg_aterm_correction, cu::DeviceMemory& d_metadata,
    cu::DeviceMemory& d_subgrid) {
  const void* parameters[] = {&time_offset,    &grid_size,
                              &subgrid_size,   &image_size,
                              &w_step,         &nr_channels,
                              &nr_stations,    d_uvw,
                              d_wavenumbers,   d_visibilities,
                              d_spheroidal,    d_aterm,
                              d_aterm_indices, d_avg_aterm_correction,
                              d_metadata,      d_subgrid};

  dim3 grid(nr_subgrids);
  dim3 block(block_gridder);
535 536
  UpdateData* data = get_update_data(get_event(), powerSensor, report,
                                     &Report::update_gridder);
537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559
  start_measurement(data);
#if ENABLE_REPEAT_KERNELS
  for (int i = 0; i < NR_REPETITIONS_GRIDDER; i++)
#endif
    executestream->launchKernel(*function_gridder, grid, block, 0, parameters);
  end_measurement(data);
}

void InstanceCUDA::launch_degridder(
    int time_offset, int nr_subgrids, int grid_size, int subgrid_size,
    float image_size, float w_step, int nr_channels, int nr_stations,
    cu::DeviceMemory& d_uvw, cu::DeviceMemory& d_wavenumbers,
    cu::DeviceMemory& d_visibilities, cu::DeviceMemory& d_spheroidal,
    cu::DeviceMemory& d_aterm, cu::DeviceMemory& d_aterm_indices,
    cu::DeviceMemory& d_metadata, cu::DeviceMemory& d_subgrid) {
  const void* parameters[] = {&time_offset,    &grid_size,   &subgrid_size,
                              &image_size,     &w_step,      &nr_channels,
                              &nr_stations,    d_uvw,        d_wavenumbers,
                              d_visibilities,  d_spheroidal, d_aterm,
                              d_aterm_indices, d_metadata,   d_subgrid};

  dim3 grid(nr_subgrids);
  dim3 block(block_degridder);
560 561
  UpdateData* data = get_update_data(get_event(), powerSensor, report,
                                     &Report::update_degridder);
562 563 564 565 566 567 568 569 570
  start_measurement(data);
#if ENABLE_REPEAT_KERNELS
  for (int i = 0; i < NR_REPETITIONS_GRIDDER; i++)
#endif
    executestream->launchKernel(*function_degridder, grid, block, 0,
                                parameters);
  end_measurement(data);
}

571
void InstanceCUDA::launch_average_beam(
Bram Veenboer's avatar
Bram Veenboer committed
572 573 574 575 576
    int nr_baselines, int nr_antennas, int nr_timesteps, int nr_channels,
    int nr_aterms, int subgrid_size, cu::DeviceMemory& d_uvw,
    cu::DeviceMemory& d_baselines, cu::DeviceMemory& d_aterms,
    cu::DeviceMemory& d_aterms_offsets, cu::DeviceMemory& d_weights,
    cu::DeviceMemory& d_average_beam) {
577 578 579 580 581 582 583 584
  const void* parameters[] = {&nr_antennas, &nr_timesteps, &nr_channels,
                              &nr_aterms,   &subgrid_size, d_uvw,
                              d_baselines,  d_aterms,      d_aterms_offsets,
                              d_weights,    d_average_beam};

  dim3 grid(nr_baselines);
  dim3 block(128);

585 586 587
  UpdateData* data = get_update_data(get_event(), powerSensor, report,
                                     &Report::update_average_beam);
  start_measurement(data);
588 589
  executestream->launchKernel(*function_average_beam, grid, block, 0,
                              parameters);
590
  end_measurement(data);
591 592
}

593 594 595 596 597 598 599 600 601 602 603 604 605
void InstanceCUDA::launch_calibrate(
    int nr_subgrids, int grid_size, int subgrid_size, float image_size,
    float w_step, int total_nr_timesteps, int nr_channels, int nr_stations,
    int nr_terms, cu::DeviceMemory& d_uvw, cu::DeviceMemory& d_wavenumbers,
    cu::DeviceMemory& d_visibilities, cu::DeviceMemory& d_weights,
    cu::DeviceMemory& d_aterm, cu::DeviceMemory& d_aterm_derivatives,
    cu::DeviceMemory& d_aterm_indices, cu::DeviceMemory& d_metadata,
    cu::DeviceMemory& d_subgrid, cu::DeviceMemory& d_sums1,
    cu::DeviceMemory& d_sums2, cu::DeviceMemory& d_lmnp,
    cu::DeviceMemory& d_hessian, cu::DeviceMemory& d_gradient,
    cu::DeviceMemory& d_residual) {
  dim3 grid(nr_subgrids);
  dim3 block(block_calibrate);
606 607
  UpdateData* data = get_update_data(get_event(), powerSensor, report,
                                     &Report::update_calibrate);
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 683 684 685 686 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 726 727 728 729 730 731 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
  start_measurement(data);

  // Get functions
  std::unique_ptr<cu::Function>& function_lmnp = functions_calibrate[0];
  std::unique_ptr<cu::Function>& function_sums = functions_calibrate[1];
  std::unique_ptr<cu::Function>& function_gradient = functions_calibrate[2];
  std::unique_ptr<cu::Function>& function_hessian = functions_calibrate[3];

  // Precompute l,m,n and phase offset
  const void* parameters_lmnp[] = {&grid_size, &subgrid_size, &image_size,
                                   &w_step,    &d_metadata,   &d_lmnp};
  executestream->launchKernel(*function_lmnp, grid, block, 0, parameters_lmnp);

  unsigned int max_nr_terms = 8;
  unsigned int current_nr_terms_y = max_nr_terms;
  for (unsigned int term_offset_y = 0; term_offset_y < (unsigned int)nr_terms;
       term_offset_y += current_nr_terms_y) {
    unsigned int last_term_y =
        min(nr_terms, term_offset_y + current_nr_terms_y);
    unsigned int current_nr_terms_y = last_term_y - term_offset_y;

    // Compute sums1
    const void* parameters_sums[] = {&subgrid_size,
                                     &image_size,
                                     &total_nr_timesteps,
                                     &nr_channels,
                                     &nr_stations,
                                     &term_offset_y,
                                     &current_nr_terms_y,
                                     &nr_terms,
                                     d_uvw,
                                     d_wavenumbers,
                                     d_aterm,
                                     d_aterm_derivatives,
                                     d_aterm_indices,
                                     d_metadata,
                                     d_subgrid,
                                     d_sums1,
                                     d_lmnp};
    executestream->launchKernel(*function_sums, grid, block, 0,
                                parameters_sums);

    // Compute gradient (diagonal)
    if (term_offset_y == 0) {
      const void* parameters_gradient[] = {&subgrid_size,
                                           &image_size,
                                           &total_nr_timesteps,
                                           &nr_channels,
                                           &nr_stations,
                                           &term_offset_y,
                                           &current_nr_terms_y,
                                           &nr_terms,
                                           d_uvw,
                                           d_wavenumbers,
                                           d_visibilities,
                                           d_weights,
                                           d_aterm,
                                           d_aterm_derivatives,
                                           d_aterm_indices,
                                           d_metadata,
                                           d_subgrid,
                                           d_sums1,
                                           d_lmnp,
                                           d_gradient,
                                           d_residual};
      executestream->launchKernel(*function_gradient, grid, block, 0,
                                  parameters_gradient);
    }

    // Compute hessian (diagonal)
    const void* parameters_hessian1[] = {&total_nr_timesteps,
                                         &nr_channels,
                                         &term_offset_y,
                                         &term_offset_y,
                                         &nr_terms,
                                         d_weights,
                                         d_aterm_indices,
                                         d_metadata,
                                         d_sums1,
                                         d_sums1,
                                         d_hessian};
    dim3 block_hessian(current_nr_terms_y, current_nr_terms_y);
    executestream->launchKernel(*function_hessian, grid, block_hessian, 0,
                                parameters_hessian1);

    unsigned int current_nr_terms_x = max_nr_terms;
    for (unsigned int term_offset_x = last_term_y;
         term_offset_x < (unsigned int)nr_terms;
         term_offset_x += current_nr_terms_x) {
      unsigned int last_term_x =
          min(nr_terms, term_offset_x + current_nr_terms_x);
      current_nr_terms_x = last_term_x - term_offset_x;

      // Compute sums2 (horizontal offset)
      const void* parameters_sums[] = {&subgrid_size,
                                       &image_size,
                                       &total_nr_timesteps,
                                       &nr_channels,
                                       &nr_stations,
                                       &term_offset_x,
                                       &current_nr_terms_x,
                                       &nr_terms,
                                       d_uvw,
                                       d_wavenumbers,
                                       d_aterm,
                                       d_aterm_derivatives,
                                       d_aterm_indices,
                                       d_metadata,
                                       d_subgrid,
                                       d_sums2,
                                       d_lmnp};
      executestream->launchKernel(*function_sums, grid, block, 0,
                                  parameters_sums);

      // Compute gradient (horizontal offset)
      if (term_offset_y == 0) {
        const void* parameters_gradient[] = {&subgrid_size,
                                             &image_size,
                                             &total_nr_timesteps,
                                             &nr_channels,
                                             &nr_stations,
                                             &term_offset_x,
                                             &current_nr_terms_x,
                                             &nr_terms,
                                             d_uvw,
                                             d_wavenumbers,
                                             d_visibilities,
                                             d_weights,
                                             d_aterm,
                                             d_aterm_derivatives,
                                             d_aterm_indices,
                                             d_metadata,
                                             d_subgrid,
                                             d_sums2,
                                             d_lmnp,
                                             d_gradient};
        executestream->launchKernel(*function_gradient, grid, block, 0,
                                    parameters_gradient);
      }

      // Compute hessian (horizontal offset)
      const void* parameters_hessian2[] = {&total_nr_timesteps,
                                           &nr_channels,
                                           &term_offset_y,
                                           &term_offset_x,
                                           &nr_terms,
                                           d_weights,
                                           d_aterm_indices,
                                           d_metadata,
                                           d_sums1,
                                           d_sums2,
                                           d_hessian};
      dim3 block_hessian(current_nr_terms_x, current_nr_terms_y);
      executestream->launchKernel(*function_hessian, grid, block_hessian, 0,
                                  parameters_hessian2);
    }
  }
  end_measurement(data);
}

void InstanceCUDA::launch_grid_fft(cu::DeviceMemory& d_data, int grid_size,
                                   DomainAtoDomainB direction) {
770
  cu::ScopedContext scc(*context);
771

772 773 774 775 776 777 778 779 780 781 782 783 784 785
  int sign =
      (direction == FourierDomainToImageDomain) ? CUFFT_INVERSE : CUFFT_FORWARD;

  // Plan FFT
  if (grid_size != m_fft_grid_size) {
    if (m_fft_plan_grid) {
      m_fft_plan_grid.reset();
    }
    m_fft_grid_size = grid_size;
    m_fft_plan_grid.reset(new cufft::C2C_2D(grid_size, grid_size));
    m_fft_plan_grid->setStream(*executestream);
  }

  // Enqueue start of measurement
786 787
  UpdateData* data = get_update_data(get_event(), powerSensor, report,
                                     &Report::update_grid_fft);
788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836
  start_measurement(data);

#if ENABLE_REPEAT_KERNELS
  for (int i = 0; i < NR_REPETITIONS_GRID_FFT; i++) {
#endif

    // Enqueue fft for every correlation
    for (unsigned i = 0; i < NR_CORRELATIONS; i++) {
      cufftComplex* data_ptr =
          reinterpret_cast<cufftComplex*>(static_cast<CUdeviceptr>(d_data));
      data_ptr += i * grid_size * grid_size;
      m_fft_plan_grid->execute(data_ptr, data_ptr, sign);
    }

#if ENABLE_REPEAT_KERNELS
  }
#endif

  // Enqueue end of measurement
  end_measurement(data);
}

void InstanceCUDA::plan_subgrid_fft(unsigned size, unsigned batch) {
#if USE_CUSTOM_FFT
  if (size == 32) {
    m_fft_subgrid_size = size;
    m_fft_subgrid_bulk = batch;
    return;
  }
#endif

  // Force plan (re-)creation if subgrid size changed
  if (size != m_fft_subgrid_size) {
    m_fft_subgrid_bulk = m_fft_subgrid_bulk_default;
    m_fft_plan_subgrid.reset();
    m_fft_subgrid_size = size;
  }

  while (m_fft_plan_subgrid == nullptr) {
    try {
      // Plan bulk fft
      unsigned stride = 1;
      unsigned dist = size * size;
      auto fft_plan = new cufft::C2C_2D(size, size, stride, dist,
                                        m_fft_subgrid_bulk * NR_CORRELATIONS);
      fft_plan->setStream(*executestream);
      m_fft_plan_subgrid.reset(fft_plan);
      auto sizeof_subgrids =
          auxiliary::sizeof_subgrids(m_fft_subgrid_bulk, m_fft_subgrid_size);
837
      d_fft_subgrid.reset(new cu::DeviceMemory(*context, sizeof_subgrids));
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
    } catch (cufft::Error& e) {
      // bulk might be too large, try again using half the bulk size
      m_fft_subgrid_bulk /= 2;
      if (m_fft_subgrid_bulk > 0) {
        std::clog << __func__ << ": reducing subgrid-fft bulk size to: "
                  << m_fft_subgrid_bulk << std::endl;
      } else {
        std::cerr << __func__ << ": could not plan subgrid-fft." << std::endl;
        throw e;
      }
    }
  }
}

void InstanceCUDA::launch_subgrid_fft(cu::DeviceMemory& d_data,
                                      unsigned int nr_subgrids,
                                      DomainAtoDomainB direction) {
  cufftComplex* data_ptr =
      reinterpret_cast<cufftComplex*>(static_cast<CUdeviceptr>(d_data));
  int sign =
      (direction == FourierDomainToImageDomain) ? CUFFT_INVERSE : CUFFT_FORWARD;

#if USE_CUSTOM_FFT
  if (fft_subgrid_size == 32) {
    const void* parameters[] = {&data_ptr, &data_ptr, &sign};
    dim3 block(128);
    dim3 grid(NR_CORRELATIONS * fft_subgrid_batch);
    executestream->launchKernel(*function_fft, grid, block, 0, parameters);
    return;
  }
#endif

870
  cu::ScopedContext scc(*context);
871

872
  // Enqueue start of measurement
873 874
  UpdateData* data = get_update_data(get_event(), powerSensor, report,
                                     &Report::update_subgrid_fft);
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
  start_measurement(data);

  // Execute bulk subgrid fft
  unsigned s = 0;
  for (; (s + m_fft_subgrid_bulk) <= nr_subgrids; s += m_fft_subgrid_bulk) {
    m_fft_plan_subgrid->execute(data_ptr, data_ptr, sign);
    data_ptr += m_fft_subgrid_size * m_fft_subgrid_size * NR_CORRELATIONS *
                m_fft_subgrid_bulk;
  }

  // Check for remainder
  unsigned int fft_subgrid_remainder = nr_subgrids % m_fft_subgrid_bulk;
  if (fft_subgrid_remainder > 0) {
    auto sizeof_subgrids =
        auxiliary::sizeof_subgrids(fft_subgrid_remainder, m_fft_subgrid_size);
    executestream->memcpyDtoDAsync(*d_fft_subgrid, (CUdeviceptr)data_ptr,
                                   sizeof_subgrids);
    cufftComplex* tmp_ptr = reinterpret_cast<cufftComplex*>(
        static_cast<CUdeviceptr>(*d_fft_subgrid));
    m_fft_plan_subgrid->execute(tmp_ptr, tmp_ptr, sign);
    executestream->memcpyDtoDAsync((CUdeviceptr)data_ptr, (CUdeviceptr)tmp_ptr,
                                   sizeof_subgrids);
  }
  executestream->synchronize();

  // Enqueue end of measurement
  end_measurement(data);
}

void InstanceCUDA::launch_grid_fft_unified(unsigned long size,
                                           unsigned int batch,
                                           Array3D<std::complex<float>>& grid,
                                           DomainAtoDomainB direction) {
908
  cu::ScopedContext scc(*context);
909

910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929
  int sign =
      (direction == FourierDomainToImageDomain) ? CUFFT_INVERSE : CUFFT_FORWARD;
  cufft::C2C_1D fft_plan_row(size, 1, 1, 1);
  cufft::C2C_1D fft_plan_col(size, size, 1, 1);

  for (unsigned i = 0; i < batch; i++) {
    // Execute 1D FFT over all columns
    for (unsigned col = 0; col < size; col++) {
      cufftComplex* ptr = (cufftComplex*)grid.data(i, col, 0);
      fft_plan_row.execute(ptr, ptr, sign);
    }

    // Execute 1D FFT over all rows
    for (unsigned row = 0; row < size; row++) {
      cufftComplex* ptr = (cufftComplex*)grid.data(i, 0, row);
      fft_plan_col.execute(ptr, ptr, sign);
    }
  }
}

Bram Veenboer's avatar
Bram Veenboer committed
930 931
void InstanceCUDA::launch_fft_shift(cu::DeviceMemory& d_data, int batch,
                                    long size, std::complex<float> scale) {
Bram Veenboer's avatar
Bram Veenboer committed
932 933
  const void* parameters[] = {&size, d_data, &scale};

Bram Veenboer's avatar
Bram Veenboer committed
934
  dim3 grid(batch, ceil(size / 2.0));
Bram Veenboer's avatar
Bram Veenboer committed
935 936
  dim3 block(128);

Bram Veenboer's avatar
Bram Veenboer committed
937 938
  UpdateData* data = get_update_data(get_event(), powerSensor, report,
                                     &Report::update_fft_shift);
Bram Veenboer's avatar
Bram Veenboer committed
939
  start_measurement(data);
Bram Veenboer's avatar
Bram Veenboer committed
940
  executestream->launchKernel(*function_fft_shift, grid, block, 0, parameters);
Bram Veenboer's avatar
Bram Veenboer committed
941 942 943
  end_measurement(data);
}

944 945 946 947 948 949 950 951 952
void InstanceCUDA::launch_adder(int nr_subgrids, long grid_size,
                                int subgrid_size, cu::DeviceMemory& d_metadata,
                                cu::DeviceMemory& d_subgrid,
                                cu::DeviceMemory& d_grid) {
  const bool enable_tiling = false;
  const void* parameters[] = {&grid_size, &subgrid_size, d_metadata,
                              d_subgrid,  d_grid,        &enable_tiling};
  dim3 grid(nr_subgrids);
  UpdateData* data =
953
      get_update_data(get_event(), powerSensor, report, &Report::update_adder);
954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972
  start_measurement(data);
#if ENABLE_REPEAT_KERNELS
  for (int i = 0; i < NR_REPETITIONS_ADDER; i++)
#endif
    executestream->launchKernel(*function_adder, grid, block_adder, 0,
                                parameters);
  end_measurement(data);
}

void InstanceCUDA::launch_adder_unified(int nr_subgrids, long grid_size,
                                        int subgrid_size,
                                        cu::DeviceMemory& d_metadata,
                                        cu::DeviceMemory& d_subgrid,
                                        void* u_grid) {
  const bool enable_tiling = true;
  const void* parameters[] = {&grid_size, &subgrid_size, d_metadata,
                              d_subgrid,  &u_grid,       &enable_tiling};
  dim3 grid(nr_subgrids);
  UpdateData* data =
973
      get_update_data(get_event(), powerSensor, report, &Report::update_adder);
974 975 976 977 978 979 980 981 982 983 984 985 986 987 988
  start_measurement(data);
  executestream->launchKernel(*function_adder, grid, block_adder, 0,
                              parameters);
  end_measurement(data);
}

void InstanceCUDA::launch_splitter(int nr_subgrids, long grid_size,
                                   int subgrid_size,
                                   cu::DeviceMemory& d_metadata,
                                   cu::DeviceMemory& d_subgrid,
                                   cu::DeviceMemory& d_grid) {
  const bool enable_tiling = false;
  const void* parameters[] = {&grid_size, &subgrid_size, d_metadata,
                              d_subgrid,  d_grid,        &enable_tiling};
  dim3 grid(nr_subgrids);
989 990
  UpdateData* data = get_update_data(get_event(), powerSensor, report,
                                     &Report::update_splitter);
991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008
  start_measurement(data);
#if ENABLE_REPEAT_KERNELS
  for (int i = 0; i < NR_REPETITIONS_ADDER; i++)
#endif
    executestream->launchKernel(*function_splitter, grid, block_splitter, 0,
                                parameters);
  end_measurement(data);
}

void InstanceCUDA::launch_splitter_unified(int nr_subgrids, long grid_size,
                                           int subgrid_size,
                                           cu::DeviceMemory& d_metadata,
                                           cu::DeviceMemory& d_subgrid,
                                           void* u_grid) {
  const bool enable_tiling = true;
  const void* parameters[] = {&grid_size, &subgrid_size, d_metadata,
                              d_subgrid,  &u_grid,       &enable_tiling};
  dim3 grid(nr_subgrids);
1009 1010
  UpdateData* data = get_update_data(get_event(), powerSensor, report,
                                     &Report::update_splitter);
1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021
  start_measurement(data);
  executestream->launchKernel(*function_splitter, grid, block_splitter, 0,
                              parameters);
  end_measurement(data);
}

void InstanceCUDA::launch_scaler(int nr_subgrids, int subgrid_size,
                                 cu::DeviceMemory& d_subgrid) {
  const void* parameters[] = {&subgrid_size, d_subgrid};
  dim3 grid(nr_subgrids);
  UpdateData* data =
1022
      get_update_data(get_event(), powerSensor, report, &Report::update_scaler);
1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080
  start_measurement(data);
  executestream->launchKernel(*function_scaler, grid, block_scaler, 0,
                              parameters);
  end_measurement(data);
}

typedef struct {
  int nr_timesteps;
  int nr_subgrids;
  Report* report;
} ReportData;

void report_job(CUstream, CUresult, void* userData) {
  ReportData* data = static_cast<ReportData*>(userData);
  int nr_timesteps = data->nr_timesteps;
  int nr_subgrids = data->nr_subgrids;
  Report* report = data->report;
  report->print(nr_timesteps, nr_subgrids);
  delete data;
}

ReportData* get_report_data(int nr_timesteps, int nr_subgrids, Report* report) {
  ReportData* data = new ReportData();
  data->nr_timesteps = nr_timesteps;
  data->nr_subgrids = nr_subgrids;
  data->report = report;
  return data;
}

void InstanceCUDA::enqueue_report(cu::Stream& stream, int nr_timesteps,
                                  int nr_subgrids) {
  ReportData* data = get_report_data(nr_timesteps, nr_subgrids, report);
  stream.addCallback((CUstreamCallback)&report_job, data);
}

void InstanceCUDA::copy_htoh(void* dst, void* src, size_t bytes) {
  char message[80];
  snprintf(message, 80, "copy_htoh(%lu)", bytes);
  cu::Marker marker(message, cu::Marker::red);
  marker.start();
  size_t batch = 1024 * 1024 * 1024;  // 1024 Mb
#pragma omp parallel for
  for (size_t i = 0; i < bytes; i += batch) {
    size_t n = i + batch < bytes ? batch : bytes - i;
    size_t src_ptr = (size_t)src + i;
    size_t dst_ptr = (size_t)dst + i;
    memcpy((void*)dst_ptr, (void*)src_ptr, n);
  }
  marker.end();
}

void InstanceCUDA::copy_dtoh(cu::Stream& stream, void* dst,
                             cu::DeviceMemory& src, size_t bytes) {
  char message[80];
  snprintf(message, 80, "copy_dtoh(%lu)", bytes);
  cu::Marker marker(message, cu::Marker::red);
  marker.start();
  size_t batch = 1024 * 1024 * 1024;  // 1024 Mb
1081
  cu::HostMemory tmp(*context, batch);
1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099
  for (size_t i = 0; i < bytes; i += batch) {
    size_t n = i + batch < bytes ? batch : bytes - i;
    size_t src_ptr = (size_t)src + i;
    size_t dst_ptr = (size_t)dst + i;
    stream.memcpyDtoHAsync((void*)tmp, (CUdeviceptr)src_ptr, n);
    stream.synchronize();
    memcpy((void*)dst_ptr, tmp, n);
  }
  marker.end();
}

void InstanceCUDA::copy_htod(cu::Stream& stream, cu::DeviceMemory& dst,
                             void* src, size_t bytes) {
  char message[80];
  snprintf(message, 80, "copy_htod(%lu)", bytes);
  cu::Marker marker(message, cu::Marker::red);
  marker.start();
  size_t batch = 1024 * 1024 * 1024;  // 1024 Mb
1100
  cu::HostMemory tmp(*context, batch);
1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118
  for (size_t i = 0; i < bytes; i += batch) {
    size_t n = i + batch < bytes ? batch : bytes - i;
    size_t src_ptr = (size_t)src + i;
    size_t dst_ptr = (size_t)dst + i;
    memcpy(tmp, (void*)src_ptr, n);
    stream.memcpyHtoDAsync((CUdeviceptr)dst_ptr, tmp, n);
    stream.synchronize();
  }
  marker.end();
}

/*
 *  Memory management per device
 *      Maintain one memory object per data structure
 */
template <typename T>
T* InstanceCUDA::reuse_memory(uint64_t size, std::unique_ptr<T>& memory) {
  if (!memory) {
1119
    memory.reset(new T(*context, size));
1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173
  } else {
    memory->resize(size);
  }
  return memory.get();
}

cu::DeviceMemory& InstanceCUDA::allocate_device_grid(size_t bytes) {
  return *reuse_memory(bytes, d_grid);
}

cu::HostMemory& InstanceCUDA::allocate_host_visibilities(size_t bytes) {
  return *reuse_memory(bytes, h_visibilities);
}

cu::HostMemory& InstanceCUDA::allocate_host_subgrids(size_t bytes) {
  return *reuse_memory(bytes, h_subgrids);
}

cu::HostMemory& InstanceCUDA::allocate_host_uvw(size_t bytes) {
  return *reuse_memory(bytes, h_uvw);
}

cu::DeviceMemory& InstanceCUDA::allocate_device_aterms(size_t bytes) {
  return *reuse_memory(bytes, d_aterms);
}

cu::DeviceMemory& InstanceCUDA::allocate_device_aterms_indices(size_t bytes) {
  return *reuse_memory(bytes, d_aterms_indices);
}

cu::DeviceMemory& InstanceCUDA::allocate_device_wavenumbers(size_t bytes) {
  return *reuse_memory(bytes, d_wavenumbers);
}

cu::DeviceMemory& InstanceCUDA::allocate_device_spheroidal(size_t bytes) {
  return *reuse_memory(bytes, d_spheroidal);
}

cu::DeviceMemory& InstanceCUDA::allocate_device_avg_aterm_correction(
    size_t bytes) {
  return *reuse_memory(bytes, d_avg_aterm_correction);
}

/*
 *  Memory management per stream
 *      Maintain multiple memory objects per structure
 *      Automatically increases the internal vectors for these objects
 */
template <typename T>
T* InstanceCUDA::reuse_memory(std::vector<std::unique_ptr<T>>& memories,
                              unsigned int id, uint64_t size) {
  T* ptr = NULL;

  if (memories.size() <= id) {
1174
    ptr = new T(*context, size);
1175

1176
    memories.emplace_back(ptr);
1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212
  } else {
    ptr = memories[id].get();
  }

  ptr->resize(size);

  return ptr;
}

cu::DeviceMemory& InstanceCUDA::allocate_device_visibilities(unsigned int id,
                                                             size_t bytes) {
  return *reuse_memory(d_visibilities_, id, bytes);
}

cu::DeviceMemory& InstanceCUDA::allocate_device_uvw(unsigned int id,
                                                    size_t bytes) {
  return *reuse_memory(d_uvw_, id, bytes);
}

cu::DeviceMemory& InstanceCUDA::allocate_device_subgrids(unsigned int id,
                                                         size_t bytes) {
  return *reuse_memory(d_subgrids_, id, bytes);
}

cu::DeviceMemory& InstanceCUDA::allocate_device_metadata(unsigned int id,
                                                         size_t bytes) {
  return *reuse_memory(d_metadata_, id, bytes);
}

/*
 * Memory management for misc device buffers
 *      Rather than storing these buffers by name,
 *      the caller gets an id that is also used to retrieve
 *      the memory from the d_misc_ vector
 */
unsigned int InstanceCUDA::allocate_device_memory(size_t bytes) {
1213
  d_misc_.emplace_back(new cu::DeviceMemory(*context, bytes));
1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234
  return d_misc_.size() - 1;
}

cu::DeviceMemory& InstanceCUDA::retrieve_device_memory(unsigned int id) {
  return *d_misc_[id];
}

/*
 * Memory management for misc page-locked host buffers
 *      Page-locking arbitrary buffers is potentially very dangerous
 *      as buffers may (partially) overlap, which will result in CUDA
 *      errors. This mechanism should only be used to register buffers
 *      that are guaranteed to be distinct and have a lifetime longer
 *      than the current InstanceCUDA object.
 */
void InstanceCUDA::register_host_memory(void* ptr, size_t bytes) {
  for (auto& memory : h_registered_) {
    if (memory->ptr() == ptr && memory->size() == bytes) {
      return;
    }
  }
1235
  h_registered_.emplace_back(new cu::RegisteredMemory(*context, ptr, bytes));
1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248 1249 1250 1251 1252 1253 1254 1255 1256 1257 1258 1259 1260 1261 1262 1263 1264 1265
}

/*
 * Host memory destructor
 */
void InstanceCUDA::free_host_memory() {
  h_visibilities.reset();
  h_uvw.reset();
  h_subgrids.reset();
  h_registered_.clear();
}

/*
 * Device memory destructor
 */
void InstanceCUDA::free_device_memory() {
  d_visibilities_.clear();
  d_uvw_.clear();
  d_metadata_.clear();
  d_subgrids_.clear();
  d_misc_.clear();
  d_aterms.reset();
  d_aterms_indices.reset();
  d_aterms_derivatives.reset();
  d_avg_aterm_correction.reset();
  d_wavenumbers.reset();
  d_spheroidal.reset();
  d_grid.reset();
}

1266 1267 1268 1269 1270
/*
 * Event destructor
 */
void InstanceCUDA::free_events() { events.clear(); }

1271 1272 1273 1274 1275 1276 1277 1278 1279 1280 1281 1282 1283 1284 1285 1286 1287 1288 1289 1290
/*
 * FFT plan destructor
 */
void InstanceCUDA::free_fft_plans() {
  m_fft_grid_size = 0;
  m_fft_plan_grid.reset();
  m_fft_subgrid_bulk = 0;
  m_fft_subgrid_size = 0;
  m_fft_plan_subgrid.reset();
  d_fft_subgrid.reset();
}

/*
 * Reset device
 */
void InstanceCUDA::reset() {
  executestream.reset();
  htodstream.reset();
  dtohstream.reset();
  context.reset(new cu::Context(*device));
1291 1292 1293
  executestream.reset(new cu::Stream(*context));
  htodstream.reset(new cu::Stream(*context));
  dtohstream.reset(new cu::Stream(*context));
1294 1295
}

1296 1297 1298 1299
/*
 * Device interface
 */
void InstanceCUDA::print_device_memory_info() const {
1300 1301 1302
#if defined(DEBUG)
  std::cout << "InstanceCUDA::" << __func__ << std::endl;
#endif
1303
  cu::ScopedContext scc(*context);
1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314
  auto memory_total =
      device->get_total_memory() / ((float)1024 * 1024 * 1024);  // GBytes
  auto memory_free =
      device->get_free_memory() / ((float)1024 * 1024 * 1024);  // GBytes
  auto memory_used = memory_total - memory_free;
  std::clog << "Device memory -> ";
  std::clog << "total: " << memory_total << " Gb, ";
  std::clog << "used: " << memory_used << " Gb, ";
  std::clog << "free: " << memory_free << " Gb" << std::endl;
}

1315
size_t InstanceCUDA::get_free_memory() const {
1316
  cu::ScopedContext scc(*context);
1317 1318 1319 1320
  return device->get_free_memory();
}

size_t InstanceCUDA::get_total_memory() const {
1321
  cu::ScopedContext scc(*context);
1322 1323 1324 1325 1326
  return device->get_total_memory();
}

template <CUdevice_attribute attribute>
int InstanceCUDA::get_attribute() const {
1327
  cu::ScopedContext scc(*context);
1328 1329 1330
  return device->get_attribute<attribute>();
}

1331 1332 1333
}  // end namespace cuda
}  // end namespace kernel
}  // end namespace idg