1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 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 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94
|
/*
//
// Copyright 2010 SRI International
//
// This file is part of the Computational Morphometry Toolkit.
//
// http://www.nitrc.org/projects/cmtk/
//
// The Computational Morphometry Toolkit is free software: you can
// redistribute it and/or modify it under the terms of the GNU General Public
// License as published by the Free Software Foundation, either version 3 of
// the License, or (at your option) any later version.
//
// The Computational Morphometry Toolkit 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 General Public License for more details.
//
// You should have received a copy of the GNU General Public License along
// with the Computational Morphometry Toolkit. If not, see
// <http://www.gnu.org/licenses/>.
//
// $Revision: 5436 $
//
// $LastChangedDate: 2018-12-10 19:01:20 -0800 (Mon, 10 Dec 2018) $
//
// $LastChangedBy: torstenrohlfing $
//
*/
#include "cmtkDeviceArrayCUDA.h"
#include <cuda_runtime_api.h>
#include <cuda_runtime.h>
#include <GPU/cmtkCUDA.h>
cmtk::DeviceArrayCUDA
::DeviceArrayCUDA( const FixedVector<3,int>& dims3 )
: m_Dims( dims3 )
{
const struct cudaChannelFormatDesc desc = cudaCreateChannelDesc<float>();
cmtkCheckCallCUDA( cudaMalloc3DArray( &(this->m_DeviceArrayPtr), &desc, make_cudaExtent( this->m_Dims[0], this->m_Dims[1], this->m_Dims[2] ) ) );
}
cmtk::DeviceArrayCUDA
::~DeviceArrayCUDA()
{
if ( this->m_DeviceArrayPtr )
cmtkCheckCallCUDA( cudaFreeArray( this->m_DeviceArrayPtr ) );
}
void
cmtk::DeviceArrayCUDA
::CopyToDevice( const float* data )
{
cudaMemcpy3DParms copyParams = {0};
copyParams.srcPtr = make_cudaPitchedPtr( (void*)data, this->m_Dims[0]*sizeof(float), this->m_Dims[0], this->m_Dims[1] );
copyParams.dstArray = this->m_DeviceArrayPtr;
copyParams.extent = make_cudaExtent( this->m_Dims[0], this->m_Dims[1], this->m_Dims[2] );
copyParams.kind = cudaMemcpyHostToDevice;
cmtkCheckCallCUDA( cudaMemcpy3D( ©Params ) );
}
void
cmtk::DeviceArrayCUDA
::CopyOnDeviceToArray( const float* data )
{
cudaMemcpy3DParms copyParams = {0};
copyParams.srcPtr = make_cudaPitchedPtr( (void*)data, this->m_Dims[0]*sizeof(float), this->m_Dims[0], this->m_Dims[1] );
copyParams.dstArray = this->m_DeviceArrayPtr;
copyParams.extent = make_cudaExtent( this->m_Dims[0], this->m_Dims[1], this->m_Dims[2] );
copyParams.kind = cudaMemcpyDeviceToDevice;
cmtkCheckCallCUDA( cudaMemcpy3D( ©Params ) );
}
void
cmtk::DeviceArrayCUDA
::CopyOnDeviceToLinear( float* data )
{
cudaMemcpy3DParms copyParams = {0};
copyParams.dstPtr = make_cudaPitchedPtr( (void*)data, this->m_Dims[0]*sizeof(float), this->m_Dims[0], this->m_Dims[1] );
copyParams.srcArray = this->m_DeviceArrayPtr;
copyParams.extent = make_cudaExtent( this->m_Dims[0], this->m_Dims[1], this->m_Dims[2] );
copyParams.kind = cudaMemcpyDeviceToDevice;
cmtkCheckCallCUDA( cudaMemcpy3D( ©Params ) );
}
|