1// FEAT3: Finite Element Analysis Toolbox, Version 3
 
    2// Copyright (C) 2010 by Stefan Turek & the FEAT group
 
    3// FEAT3 is released under the GNU General Public License version 3,
 
    4// see the file 'copyright.txt' in the top level directory for details.
 
    7#include <kernel/base_header.hpp>
 
    8#include <kernel/lafem/arch/component_copy.hpp>
 
    9#include <kernel/util/exception.hpp>
 
   10#include <kernel/util/memory_pool.hpp>
 
   11#include <kernel/util/cuda_util.hpp>
 
   12#include <kernel/util/half.hpp>
 
   18using namespace FEAT::LAFEM;
 
   19using namespace FEAT::LAFEM::Arch;
 
   21void ComponentCopy::value_cuda(float * r, const float * const x, const int stride, const int block, const Index size)
 
   23  cublasCopyEx(Util::Intern::cublas_handle, int(size), x, CUDA_R_32F, 1, &r[block], CUDA_R_32F, stride);
 
   24  cudaDeviceSynchronize();
 
   26  cudaError_t last_error(cudaGetLastError());
 
   27  if (cudaSuccess != last_error)
 
   28    throw InternalError(__func__, __FILE__, __LINE__, "CUDA error occurred in execution!\n" + stringify(cudaGetErrorString(last_error)));
 
   32void ComponentCopy::value_to_cuda(const float * const r, float * x, const int stride, const int block, const Index size)
 
   34  cublasCopyEx(Util::Intern::cublas_handle, int(size), &r[block], CUDA_R_32F, stride, x, CUDA_R_32F, 1);
 
   35  cudaDeviceSynchronize();
 
   37  cudaError_t last_error(cudaGetLastError());
 
   38  if (cudaSuccess != last_error)
 
   39    throw InternalError(__func__, __FILE__, __LINE__, "CUDA error occurred in execution!\n" + stringify(cudaGetErrorString(last_error)));
 
   43void ComponentCopy::value_cuda(double * r, const double * const x, const int stride, const int block, const Index size)
 
   45  cublasCopyEx(Util::Intern::cublas_handle, int(size), x, CUDA_R_64F, 1, &r[block], CUDA_R_64F, stride);
 
   46  cudaDeviceSynchronize();
 
   48  cudaError_t last_error(cudaGetLastError());
 
   49  if (cudaSuccess != last_error)
 
   50    throw InternalError(__func__, __FILE__, __LINE__, "CUDA error occurred in execution!\n" + stringify(cudaGetErrorString(last_error)));
 
   54void ComponentCopy::value_to_cuda(const double * const r, double * x, const int stride, const int block, const Index size)
 
   56  cublasCopyEx(Util::Intern::cublas_handle, int(size), &r[block], CUDA_R_64F, stride, x, CUDA_R_64F, 1);
 
   57  cudaDeviceSynchronize();
 
   59  cudaError_t last_error(cudaGetLastError());
 
   60  if (cudaSuccess != last_error)
 
   61    throw InternalError(__func__, __FILE__, __LINE__, "CUDA error occurred in execution!\n" + stringify(cudaGetErrorString(last_error)));
 
   65#ifdef FEAT_HAVE_HALFMATH
 
   66void ComponentCopy::value_cuda(Half * r, const Half * const x, const int stride, const int block, const Index size)
 
   68  cublasCopyEx(Util::Intern::cublas_handle, int(size), x, CUDA_R_16F, 1, &r[block], CUDA_R_16F, stride);
 
   69  cudaDeviceSynchronize();
 
   71  cudaError_t last_error(cudaGetLastError());
 
   72  if (cudaSuccess != last_error)
 
   73    throw InternalError(__func__, __FILE__, __LINE__, "CUDA error occurred in execution!\n" + stringify(cudaGetErrorString(last_error)));
 
   76void ComponentCopy::value_to_cuda(const Half * const r, Half * x, const int stride, const int block, const Index size)
 
   78  cublasCopyEx(Util::Intern::cublas_handle, int(size), &r[block], CUDA_R_16F, stride, x, CUDA_R_16F, 1);
 
   79  cudaDeviceSynchronize();
 
   81  cudaError_t last_error(cudaGetLastError());
 
   82  if (cudaSuccess != last_error)
 
   83    throw InternalError(__func__, __FILE__, __LINE__, "CUDA error occurred in execution!\n" + stringify(cudaGetErrorString(last_error)));