Add exp, mod, acos, acosd, conj, norm and unittest.
This commit is contained in:
@@ -1,14 +1,18 @@
|
||||
#include "CudaMatrix.h"
|
||||
#include "AuroraDefs.h"
|
||||
#include "Function1D.cuh"
|
||||
#include "Function1D.h"
|
||||
#include "Matrix.h"
|
||||
|
||||
#include <cmath>
|
||||
#include <cstdlib>
|
||||
#include <thrust/device_vector.h>
|
||||
#include <thrust/transform.h>
|
||||
#include <thrust/iterator/constant_iterator.h>
|
||||
#include <thrust/iterator/counting_iterator.h>
|
||||
#include <thrust/complex.h>
|
||||
#include <cuda_runtime.h>
|
||||
#include <cusolverDn.h>
|
||||
|
||||
using namespace Aurora;
|
||||
using namespace thrust::placeholders;
|
||||
@@ -701,3 +705,251 @@ CudaMatrix Aurora::log(const CudaMatrix& aMatrix, int aBaseNum)
|
||||
cudaDeviceSynchronize();
|
||||
return Aurora::CudaMatrix::fromRawData(data, aMatrix.getDimSize(0), aMatrix.getDimSize(1), aMatrix.getDimSize(2), aMatrix.getValueType());
|
||||
}
|
||||
|
||||
__global__ void expKernel(float* aInputData, float* aOutput, unsigned int aInputSize, bool aIsComplex)
|
||||
{
|
||||
unsigned int idx = blockIdx.x * blockDim.x + threadIdx.x;
|
||||
if (idx < aInputSize)
|
||||
{
|
||||
if(aIsComplex)
|
||||
{
|
||||
unsigned int index = 2 * idx;
|
||||
float expReal = expf(aInputData[index]);
|
||||
aOutput[index] = expReal * cosf(aInputData[index + 1]);
|
||||
aOutput[index + 1] = expReal * sinf(aInputData[index + 1]);
|
||||
}
|
||||
else
|
||||
{
|
||||
aOutput[idx] = expf(aInputData[idx]);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
CudaMatrix Aurora::exp(const CudaMatrix& aMatrix)
|
||||
{
|
||||
size_t size = aMatrix.getDataSize();
|
||||
float* data = nullptr;
|
||||
cudaMalloc((void**)&data, sizeof(float) * size * aMatrix.getValueType());
|
||||
int blocksPerGrid = (size + THREADS_PER_BLOCK - 1) / THREADS_PER_BLOCK;
|
||||
expKernel<<<blocksPerGrid, THREADS_PER_BLOCK>>>(aMatrix.getData(), data, size, aMatrix.isComplex());
|
||||
cudaDeviceSynchronize();
|
||||
return Aurora::CudaMatrix::fromRawData(data, aMatrix.getDimSize(0), aMatrix.getDimSize(1), aMatrix.getDimSize(2), aMatrix.getValueType());
|
||||
}
|
||||
|
||||
__global__ void modKernel(float* aInputData, float* aOutput, unsigned int aInputSize, float aModValue)
|
||||
{
|
||||
unsigned int idx = blockIdx.x * blockDim.x + threadIdx.x;
|
||||
if (idx < aInputSize)
|
||||
{
|
||||
aOutput[idx] = fmodf(aInputData[idx], aModValue);
|
||||
}
|
||||
}
|
||||
|
||||
CudaMatrix Aurora::mod(const CudaMatrix& aMatrix, float aValue)
|
||||
{
|
||||
if(aMatrix.isComplex())
|
||||
{
|
||||
std::cerr<<"mod not support complex"<<std::endl;
|
||||
return CudaMatrix();
|
||||
}
|
||||
|
||||
size_t size = aMatrix.getDataSize();
|
||||
float* data = nullptr;
|
||||
cudaMalloc((void**)&data, sizeof(float) * size);
|
||||
int blocksPerGrid = (size + THREADS_PER_BLOCK - 1) / THREADS_PER_BLOCK;
|
||||
modKernel<<<blocksPerGrid, THREADS_PER_BLOCK>>>(aMatrix.getData(), data, size, aValue);
|
||||
cudaDeviceSynchronize();
|
||||
return Aurora::CudaMatrix::fromRawData(data, aMatrix.getDimSize(0), aMatrix.getDimSize(1), aMatrix.getDimSize(2));
|
||||
}
|
||||
|
||||
__global__ void acosKernel(float* aInputData, float* aOutput, unsigned int aInputSize)
|
||||
{
|
||||
unsigned int idx = blockIdx.x * blockDim.x + threadIdx.x;
|
||||
if (idx < aInputSize)
|
||||
{
|
||||
aOutput[idx] = acosf(aInputData[idx]);
|
||||
}
|
||||
}
|
||||
|
||||
CudaMatrix Aurora::acos(const CudaMatrix& aMatrix)
|
||||
{
|
||||
if(aMatrix.isComplex())
|
||||
{
|
||||
std::cerr<<"acos not support complex"<<std::endl;
|
||||
return CudaMatrix();
|
||||
}
|
||||
|
||||
size_t size = aMatrix.getDataSize();
|
||||
float* data = nullptr;
|
||||
cudaMalloc((void**)&data, sizeof(float) * size);
|
||||
int blocksPerGrid = (size + THREADS_PER_BLOCK - 1) / THREADS_PER_BLOCK;
|
||||
acosKernel<<<blocksPerGrid, THREADS_PER_BLOCK>>>(aMatrix.getData(), data, size);
|
||||
cudaDeviceSynchronize();
|
||||
return Aurora::CudaMatrix::fromRawData(data, aMatrix.getDimSize(0), aMatrix.getDimSize(1), aMatrix.getDimSize(2));
|
||||
}
|
||||
|
||||
__global__ void acosdKernel(float* aInputData, float* aOutput, unsigned int aInputSize)
|
||||
{
|
||||
unsigned int idx = blockIdx.x * blockDim.x + threadIdx.x;
|
||||
if (idx < aInputSize)
|
||||
{
|
||||
aOutput[idx] = acosf(aInputData[idx]) * 180 / PI;
|
||||
}
|
||||
}
|
||||
|
||||
CudaMatrix Aurora::acosd(const CudaMatrix& aMatrix)
|
||||
{
|
||||
if(aMatrix.isComplex())
|
||||
{
|
||||
std::cerr<<"acos not support complex"<<std::endl;
|
||||
return CudaMatrix();
|
||||
}
|
||||
|
||||
size_t size = aMatrix.getDataSize();
|
||||
float* data = nullptr;
|
||||
cudaMalloc((void**)&data, sizeof(float) * size);
|
||||
int blocksPerGrid = (size + THREADS_PER_BLOCK - 1) / THREADS_PER_BLOCK;
|
||||
acosdKernel<<<blocksPerGrid, THREADS_PER_BLOCK>>>(aMatrix.getData(), data, size);
|
||||
cudaDeviceSynchronize();
|
||||
return Aurora::CudaMatrix::fromRawData(data, aMatrix.getDimSize(0), aMatrix.getDimSize(1), aMatrix.getDimSize(2));
|
||||
}
|
||||
|
||||
__global__ void conjKernel(float* aInputData, float* aOutput, unsigned int aInputSize)
|
||||
{
|
||||
unsigned int idx = blockIdx.x * blockDim.x + threadIdx.x;
|
||||
if (idx < aInputSize)
|
||||
{
|
||||
unsigned int index = idx * 2;
|
||||
aOutput[index] = aInputData[index];
|
||||
aOutput[index + 1] = -aInputData[index + 1];
|
||||
}
|
||||
}
|
||||
|
||||
CudaMatrix Aurora::conj(const CudaMatrix& aMatrix)
|
||||
{
|
||||
if(!aMatrix.isComplex())
|
||||
{
|
||||
return CudaMatrix::copyFromRawData(aMatrix.getData(),aMatrix.getDimSize(0),aMatrix.getDimSize(1),aMatrix.getDimSize(2));
|
||||
}
|
||||
size_t size = aMatrix.getDataSize();
|
||||
float* data = nullptr;
|
||||
cudaMalloc((void**)&data, sizeof(float) * size * aMatrix.getValueType());
|
||||
int blocksPerGrid = (size + THREADS_PER_BLOCK - 1) / THREADS_PER_BLOCK;
|
||||
conjKernel<<<blocksPerGrid, THREADS_PER_BLOCK>>>(aMatrix.getData(), data, size);
|
||||
cudaDeviceSynchronize();
|
||||
return Aurora::CudaMatrix::fromRawData(data, aMatrix.getDimSize(0), aMatrix.getDimSize(1), aMatrix.getDimSize(2), aMatrix.getValueType());
|
||||
}
|
||||
|
||||
|
||||
float Aurora::norm(const CudaMatrix& aMatrix, NormMethod aNormMethod)
|
||||
{
|
||||
float resultValue = 0;
|
||||
if(aMatrix.getDims() > 2)
|
||||
{
|
||||
std::cerr<<"norm not support 3d matrix"<<std::endl;
|
||||
return 0;
|
||||
}
|
||||
//for 1 dims
|
||||
if(aMatrix.getDimSize(0) == 1 || aMatrix.getDimSize(1) == 1 )
|
||||
{
|
||||
if(aNormMethod == Aurora::Norm1)
|
||||
{
|
||||
CudaMatrix result = abs(aMatrix);
|
||||
resultValue = thrust::reduce(thrust::device, result.getData(), result.getData() + result.getDataSize(), 0.0, thrust::plus<float>());
|
||||
return resultValue;
|
||||
}
|
||||
else
|
||||
{
|
||||
CudaMatrix result = aMatrix.deepCopy();
|
||||
thrust::transform(thrust::device, result.getData(), result.getData() + result.getDataSize() * result.getValueType(), result.getData(), thrust::square<float>());
|
||||
resultValue = thrust::reduce(thrust::device, result.getData(), result.getData() + result.getDataSize() * result.getValueType(), 0.0, thrust::plus<float>());
|
||||
return std::sqrt(resultValue);
|
||||
}
|
||||
}
|
||||
//for 2 dims
|
||||
if(aNormMethod == Aurora::NormF)
|
||||
{
|
||||
CudaMatrix result = aMatrix.deepCopy();
|
||||
thrust::transform(thrust::device, result.getData(), result.getData() + result.getDataSize() * result.getValueType(), result.getData(), thrust::square<float>());
|
||||
resultValue = thrust::reduce(thrust::device, result.getData(), result.getData() + result.getDataSize() * result.getValueType(), 0.0, thrust::plus<float>());
|
||||
return std::sqrt(resultValue);
|
||||
}
|
||||
else if(aNormMethod == Aurora::Norm1)
|
||||
{
|
||||
for(int i=0; i<aMatrix.getDimSize(1); ++i)
|
||||
{
|
||||
CudaMatrix result = abs(aMatrix.block(1, i, i));
|
||||
float tempValue = thrust::reduce(thrust::device, result.getData(), result.getData() + result.getDataSize(), 0.0, thrust::plus<float>());
|
||||
if(resultValue < tempValue)
|
||||
{
|
||||
resultValue = tempValue;
|
||||
}
|
||||
}
|
||||
return resultValue;
|
||||
}
|
||||
else
|
||||
{
|
||||
if(aMatrix.isComplex())
|
||||
{
|
||||
std::cerr<<"norm2 not support 2d complex matrix"<<std::endl;
|
||||
return 0;
|
||||
}
|
||||
cusolverDnHandle_t cusolverH = NULL;
|
||||
cudaStream_t stream = NULL;
|
||||
float* d_S = NULL;
|
||||
float* d_U = NULL;
|
||||
float* d_VT = NULL;
|
||||
int* devInfo = NULL;
|
||||
float* d_work = NULL;
|
||||
int lwork = 0;
|
||||
const int m = aMatrix.getDimSize(0);
|
||||
const int n = aMatrix.getDimSize(1);
|
||||
const int lda = m;
|
||||
const int ldu = m;
|
||||
const int ldvt = n;
|
||||
|
||||
cusolverDnCreate(&cusolverH);
|
||||
cudaStreamCreateWithFlags(&stream, cudaStreamNonBlocking);
|
||||
cusolverDnSetStream(cusolverH, stream);
|
||||
cudaMalloc((void**)&d_S, sizeof(float)*n);
|
||||
cudaMalloc((void**)&d_U, sizeof(float)*ldu*m);
|
||||
cudaMalloc((void**)&d_VT, sizeof(float)*ldvt*n);
|
||||
cudaMalloc((void**)&devInfo, sizeof(int));
|
||||
|
||||
cusolverDnSgesvd_bufferSize(cusolverH, m, n, &lwork);
|
||||
|
||||
cudaMalloc((void**)&d_work, sizeof(float)*lwork);
|
||||
auto matrix = aMatrix.deepCopy();
|
||||
cusolverDnSgesvd(cusolverH, 'A', 'A', m, n, matrix.getData(), lda, d_S,d_U, ldu, d_VT, ldvt, d_work, lwork, NULL, devInfo);
|
||||
|
||||
int devInfo_h = 0;
|
||||
cudaMemcpy(&devInfo_h, devInfo, sizeof(int), cudaMemcpyDeviceToHost);
|
||||
|
||||
if (devInfo_h != 0)
|
||||
{
|
||||
printf("Unsuccessful SVD execution\n");
|
||||
printf("Error code: %d\n", devInfo_h);
|
||||
}
|
||||
|
||||
float S[n] = {0};
|
||||
cudaMemcpy(S, d_S, sizeof(float)*n, cudaMemcpyDeviceToHost);
|
||||
|
||||
float resultValue = S[0];
|
||||
for (int i = 1; i < n; i++)
|
||||
{
|
||||
if (S[i] > resultValue)
|
||||
{
|
||||
resultValue = S[i];
|
||||
}
|
||||
}
|
||||
|
||||
if (d_S ) cudaFree(d_S);
|
||||
if (d_U ) cudaFree(d_U);
|
||||
if (d_VT) cudaFree(d_VT);
|
||||
if (devInfo) cudaFree(devInfo);
|
||||
if (d_work) cudaFree(d_work);
|
||||
if (cusolverH) cusolverDnDestroy(cusolverH);
|
||||
if (stream) cudaStreamDestroy(stream);
|
||||
return resultValue;
|
||||
}
|
||||
}
|
||||
|
||||
@@ -51,7 +51,17 @@ namespace Aurora
|
||||
|
||||
CudaMatrix log(const CudaMatrix& aMatrix, int aBaseNum = -1);
|
||||
|
||||
CudaMatrix exp(const CudaMatrix& aMatrix);
|
||||
|
||||
CudaMatrix mod(const CudaMatrix& aMatrix, float aValue);
|
||||
|
||||
CudaMatrix acos(const CudaMatrix& aMatrix);
|
||||
|
||||
CudaMatrix acosd(const CudaMatrix& aMatrix);
|
||||
|
||||
CudaMatrix conj(const CudaMatrix& aMatrix);
|
||||
|
||||
float norm(const CudaMatrix& aMatrix, NormMethod aNormMethod);
|
||||
|
||||
// ------compareSet----------------------------------------------------
|
||||
|
||||
|
||||
Reference in New Issue
Block a user