gpufilter
GPU-Efficient Recursive Filtering and Summed-Area Tables
|
00001 00009 #ifndef DVECTOR_H 00010 #define DVECTOR_H 00011 00012 //== INCLUDES ================================================================= 00013 00014 #include <vector> 00015 00016 #include <alloc.h> 00017 00018 //== NAMESPACES =============================================================== 00019 00020 namespace gpufilter { 00021 00022 //== CLASS DEFINITION ========================================================= 00023 00033 template< class T > 00034 class dvector { 00035 00036 public: 00037 00042 explicit dvector( const std::vector<T>& that ) : m_size(0), m_capacity(0), m_data(0) { 00043 *this = that; 00044 } 00045 00051 dvector( const T *data, 00052 const size_t& size ) : m_size(0), m_capacity(0), m_data(0) { 00053 copy_from( data, size ); 00054 } 00055 00064 dvector( const T *data, 00065 const size_t& w_data, 00066 const size_t& h_data, 00067 const size_t& w, 00068 const size_t& h ) : m_size(0), m_capacity(0), m_data(0) { 00069 copy_from( data, w_data, h_data, w, h ); 00070 } 00071 00076 dvector( const dvector<T>& that ) : m_size(0), m_capacity(0), m_data(0) { 00077 *this = that; 00078 } 00079 00084 dvector( const size_t& size = 0 ) : m_size(0), m_capacity(0), m_data(0) { 00085 resize(size); 00086 } 00087 00091 ~dvector() { 00092 cuda_delete(m_data); 00093 m_data = 0; 00094 m_capacity = 0; 00095 m_size = 0; 00096 } 00097 00102 void resize( const size_t& size ) { 00103 if( size > m_capacity ) 00104 { 00105 cuda_delete(m_data); 00106 m_data = 0; 00107 m_capacity = 0; 00108 m_size = 0; 00109 00110 m_data = cuda_new<T>(size); 00111 m_capacity = size; 00112 m_size = size; 00113 } 00114 else 00115 m_size = size; 00116 } 00117 00121 void clear() { 00122 m_size = 0; 00123 } 00124 00130 T operator [] ( const int& idx ) const { 00131 T value; 00132 cudaMemcpy(&value, data()+idx, sizeof(T), cudaMemcpyDeviceToHost); 00133 return value; 00134 } 00135 00141 dvector& operator = ( const dvector<T>& that ) { 00142 resize(that.size()); 00143 cudaMemcpy(data(), that.data(), size()*sizeof(T), cudaMemcpyDeviceToDevice); 00144 cuda_error("Error during memcpy from device to device"); 00145 return *this; 00146 } 00147 00153 dvector& operator = ( const std::vector<T>& that ) { 00154 resize(that.size()); 00155 cudaMemcpy(data(), &that[0], size()*sizeof(T), cudaMemcpyHostToDevice); 00156 cuda_error("Error during memcpy from host to device"); 00157 return *this; 00158 } 00159 00165 void copy_to( T *data, 00166 const size_t& s ) const { 00167 cudaMemcpy(data, this->data(), std::min(size(),s)*sizeof(T), cudaMemcpyDeviceToHost); 00168 cuda_error("Error during memcpy from device to host"); 00169 } 00170 00179 void copy_to( T *data, 00180 const size_t& w, 00181 const size_t& h, 00182 const size_t& w_data, 00183 const size_t& h_data ) const { 00184 cudaMemcpy2D(data, w_data*sizeof(T), this->data(), w*sizeof(T), w_data*sizeof(T), h_data, cudaMemcpyDeviceToHost); 00185 cuda_error("Error during memcpy2D from device to host"); 00186 } 00187 00193 void copy_from( const T *data, 00194 const size_t& size ) { 00195 resize(size); 00196 cudaMemcpy(this->data(), const_cast<T *>(data), 00197 size*sizeof(T), cudaMemcpyHostToDevice); 00198 cuda_error("Error during memcpy from host to device"); 00199 } 00200 00209 void copy_from( const T *data, 00210 const size_t& w_data, 00211 const size_t& h_data, 00212 const size_t& w, 00213 const size_t& h ) { 00214 resize(w*h); 00215 cudaMemcpy2D(this->data(), w*sizeof(T), 00216 const_cast<T *>(data), w_data*sizeof(T), w_data*sizeof(T), 00217 h_data, cudaMemcpyHostToDevice); 00218 cuda_error("Error during memcpy2D from host to device"); 00219 } 00220 00224 void fill_zero() { 00225 cudaMemset(m_data, 0, m_size*sizeof(T)); 00226 } 00227 00232 bool empty() const { return size()==0; } 00233 00238 size_t size() const { return m_size; } 00239 00244 T *data() { return m_data; } 00245 00250 const T *data() const { return m_data; } 00251 00256 T back() const { return operator[](size()-1); } 00257 00262 operator T* () { return data(); } 00263 00268 operator const T* () const { return data(); } 00269 00275 friend void swap( dvector<T>& a, 00276 dvector<T>& b ) { 00277 std::swap(a.m_data, b.m_data); 00278 std::swap(a.m_size, b.m_size); 00279 std::swap(a.m_capacity, b.m_capacity); 00280 } 00281 00282 private: 00283 00284 T *m_data; 00285 size_t m_size; 00286 size_t m_capacity; 00287 00288 }; 00289 00290 //== IMPLEMENTATION =========================================================== 00291 00303 template< class T > 00304 std::vector<T> to_cpu( const T *d_vec, 00305 unsigned len ) { 00306 std::vector<T> out; 00307 out.resize(len); 00308 00309 cudaMemcpy(&out[0], d_vec, len*sizeof(T), cudaMemcpyDeviceToHost); 00310 cuda_error("Error during memcpy from device to host"); 00311 00312 return out; 00313 } 00314 00323 template< class T > 00324 std::vector<T> to_cpu( const dvector<T>& v ) { 00325 return to_cpu(v.data(), v.size()); 00326 } 00327 00328 //============================================================================= 00329 } // namespace gpufilter 00330 //============================================================================= 00331 #endif // DVECTOR_H 00332 //============================================================================= 00333 //vi: ai sw=4 ts=4