![]() |
kspaceFirstOrder3D-CUDA
1.1
The CUDA/C++ implementation of the k-wave toolbox for the time-domain simulation of acoustic wave fields in 3D
|
The implementation file of cuda kernels used for data sampling (output streams). More...
#include <cuda.h>
#include <cuda_runtime.h>
#include <OutputHDF5Streams/BaseOutputHDF5Stream.h>
#include <OutputHDF5Streams/OutputStreamsCUDAKernels.cuh>
#include <Parameters/Parameters.h>
#include <Logger/Logger.h>
#include <Utils/CUDAUtils.cuh>
Go to the source code of this file.
Functions | |
int | GetSamplerBlockSize () |
int | GetSamplerGridSize () |
template<TBaseOutputHDF5Stream::TReduceOperator reduceOp> | |
__global__ void | CUDASampleIndex (float *samplingBuffer, const float *sourceData, const size_t *sensorData, const size_t nSamples) |
template void | OutputStreamsCUDAKernels::SampleIndex< TBaseOutputHDF5Stream::TReduceOperator::NONE > (float *samplingBuffer, const float *sourceData, const size_t *sensorData, const size_t nSamples) |
template void | OutputStreamsCUDAKernels::SampleIndex< TBaseOutputHDF5Stream::TReduceOperator::RMS > (float *samplingBuffer, const float *sourceData, const size_t *sensorData, const size_t nSamples) |
template void | OutputStreamsCUDAKernels::SampleIndex< TBaseOutputHDF5Stream::TReduceOperator::MAX > (float *samplingBuffer, const float *sourceData, const size_t *sensorData, const size_t nSamples) |
template void | OutputStreamsCUDAKernels::SampleIndex< TBaseOutputHDF5Stream::TReduceOperator::MIN > (float *samplingBuffer, const float *sourceData, const size_t *sensorData, const size_t nSamples) |
__device__ size_t | TransformCoordinates (const size_t cuboidIdx, const dim3 &topLeftCorner, const dim3 &bottomRightCorner, const dim3 &matrixSize) |
template<TBaseOutputHDF5Stream::TReduceOperator reduceOp> | |
__global__ void | CUDASampleCuboid (float *samplingBuffer, const float *sourceData, const dim3 topLeftCorner, const dim3 bottomRightCorner, const dim3 matrixSize, const size_t nSamples) |
template void | OutputStreamsCUDAKernels::SampleCuboid< TBaseOutputHDF5Stream::TReduceOperator::NONE > (float *samplingBuffer, const float *sourceData, const dim3 topLeftCorner, const dim3 bottomRightCorner, const dim3 matrixSize, const size_t nSamples) |
template void | OutputStreamsCUDAKernels::SampleCuboid< TBaseOutputHDF5Stream::TReduceOperator::RMS > (float *samplingBuffer, const float *sourceData, const dim3 topLeftCorner, const dim3 bottomRightCorner, const dim3 matrixSize, const size_t nSamples) |
template void | OutputStreamsCUDAKernels::SampleCuboid< TBaseOutputHDF5Stream::TReduceOperator::MAX > (float *samplingBuffer, const float *sourceData, const dim3 topLeftCorner, const dim3 bottomRightCorner, const dim3 matrixSize, const size_t nSamples) |
template void | OutputStreamsCUDAKernels::SampleCuboid< TBaseOutputHDF5Stream::TReduceOperator::MIN > (float *samplingBuffer, const float *sourceData, const dim3 topLeftCorner, const dim3 bottomRightCorner, const dim3 matrixSize, const size_t nSamples) |
template<TBaseOutputHDF5Stream::TReduceOperator reduceOp> | |
__global__ void | CUDASampleAll (float *samplingBuffer, const float *sourceData, const size_t nSamples) |
template void | OutputStreamsCUDAKernels::SampleAll< TBaseOutputHDF5Stream::TReduceOperator::RMS > (float *samplingBuffer, const float *sourceData, const size_t nSamples) |
template void | OutputStreamsCUDAKernels::SampleAll< TBaseOutputHDF5Stream::TReduceOperator::MAX > (float *samplingBuffer, const float *sourceData, const size_t nSamples) |
template void | OutputStreamsCUDAKernels::SampleAll< TBaseOutputHDF5Stream::TReduceOperator::MIN > (float *samplingBuffer, const float *sourceData, const size_t nSamples) |
__global__ void | CUDAPostProcessingRMS (float *samplingBuffer, const float scalingCoeff, const size_t nSamples) |
This file is part of the C++ extension of the k-Wave Toolbox (http://www.k-wave.org).
Copyright (C) 2016 Jiri Jaros and Bradley Treeby.
This file is part of the k-Wave. k-Wave is free software: you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.
k-Wave is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with k-Wave. If not, see http://www.gnu.org/licenses/.
Definition in file OutputStreamsCUDAKernels.cu.
__global__ void CUDAPostProcessingRMS | ( | float * | samplingBuffer, |
const float | scalingCoeff, | ||
const size_t | nSamples | ||
) |
CUDA kernel to apply post-processing for RMS
[in,out] | samplingBuffer | - Buffer to apply post-processing on |
[in] | scalingCoeff | - Scaling coeficinet for RMS |
[in] | nSamples | - Number of elements |
Definition at line 436 of file OutputStreamsCUDAKernels.cu.
__global__ void CUDASampleAll | ( | float * | samplingBuffer, |
const float * | sourceData, | ||
const size_t | nSamples | ||
) |
CUDA kernel to sample and aggregate the source matrix on the whole domain and apply a reduce operator.
[in,out] | samplingBuffer | - Buffer to sample data in |
[in] | sourceData | - Source matrix |
[in] | nSamples | - Number of sampled points |
Definition at line 355 of file OutputStreamsCUDAKernels.cu.
__global__ void CUDASampleCuboid | ( | float * | samplingBuffer, |
const float * | sourceData, | ||
const dim3 | topLeftCorner, | ||
const dim3 | bottomRightCorner, | ||
const dim3 | matrixSize, | ||
const size_t | nSamples | ||
) |
CUDA kernel to sample data inside one cuboid, operation is selected by a template parameter.
[out] | samplingBuffer | - Buffer to sample data in |
[in] | sourceData | - Source matrix |
[in] | topLeftCorner | - Top left corner of the cuboid |
[in] | bottomRightCorner | - Bottom right corner of the cuboid |
[in] | matrixSize | - Dimension sizes of the matrix being sampled |
[in] | nSamples | - Number of grid points inside the cuboid |
Definition at line 236 of file OutputStreamsCUDAKernels.cu.
__global__ void CUDASampleIndex | ( | float * | samplingBuffer, |
const float * | sourceData, | ||
const size_t * | sensorData, | ||
const size_t | nSamples | ||
) |
CUDA kernel to sample data based on index sensor mask. The operator is given by the template parameter.
[out] | samplingBuffer | - Buffer to sample data in |
[in] | sourceData | - Source matrix |
[in] | sensorData | - Sensor mask |
[in] | nSamples | - Number of sampled points |
Definition at line 92 of file OutputStreamsCUDAKernels.cu.
int GetSamplerBlockSize | ( | ) |
Get Sampler CUDA Block size.
Definition at line 59 of file OutputStreamsCUDAKernels.cu.
int GetSamplerGridSize | ( | ) |
Get sampler CUDA grid size.
Definition at line 71 of file OutputStreamsCUDAKernels.cu.
|
inline |
Transform 3D coordinates within the cuboid into 1D coordinates within the matrix being sampled.
[in] | cuboidIdx | - Cuboid index |
[in] | topLeftCorner | - Top left corner |
[in] | bottomRightCorner | - Bottom right corner |
[in] | matrixSize | - Size of the matrix being sampled |
Definition at line 195 of file OutputStreamsCUDAKernels.cu.