Xmipp  v3.23.11-Nereus
Classes | Macros | Functions
cuda_xmipp_utils.h File Reference
#include <stdio.h>
#include <complex>
Include dependency graph for cuda_xmipp_utils.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Classes

class  mycufftHandle
 
class  myStreamHandle
 
class  XmippDim3
 
class  TransformMatrix< T >
 
class  GpuMultidimArrayAtGpu< T >
 

Macros

#define CONVERT2DIM3(d)   (dim3((d).x,(d).y,(d).z))
 

Functions

void mycufftDestroy (void *ptr)
 
void myStreamDestroy (void *ptr)
 
void myStreamCreate (myStreamHandle &myStream)
 
void gpuMalloc (void **d_data, size_t Nbytes)
 
void gpuFree (void *d_data)
 
void cpuMalloc (void **h_data, size_t Nbytes)
 
void cpuFree (void *h_data)
 
void initializeIdentity (float *d_data, float *h_data, int Ndim, myStreamHandle &myStream)
 
void setTranslationMatrix (float *d_data, float *posX, float *posY, int Ndim, myStreamHandle &myStream)
 
void setRotationMatrix (float *d_data, float *ang, int Ndim, myStreamHandle &myStream)
 
void gpuCopyFromGPUToGPUStream (void *d_dataFrom, void *d_dataTo, size_t Nbytes, myStreamHandle &myStream)
 
void gpuCopyFromCPUToGPUStream (void *data, void *d_data, size_t Nbytes, myStreamHandle &myStream)
 
void gpuCopyFromGPUToCPUStream (void *d_data, void *data, size_t Nbytes, myStreamHandle &myStream)
 
void gpuCopyFromGPUToGPU (void *d_dataFrom, void *d_dataTo, size_t Nbytes)
 
void gpuCopyFromCPUToGPU (void *data, void *d_data, size_t Nbytes)
 
void gpuCopyFromGPUToCPU (void *d_data, void *data, size_t Nbytes)
 
int gridFromBlock (int tasks, int Nthreads)
 
template<typename T >
T * loadToGPU (const T *data, size_t items)
 
void cuda_check_gpu_memory (float *data)
 
void cuda_check_gpu_properties (int *maxGridSize)
 

Function Documentation

◆ cpuFree()

void cpuFree ( void *  h_data)

Definition at line 211 of file cuda_xmipp_utils.cpp.

212 {
213  gpuErrchk(cudaFreeHost(h_data));
214 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31

◆ cpuMalloc()

void cpuMalloc ( void **  h_data,
size_t  Nbytes 
)

Definition at line 206 of file cuda_xmipp_utils.cpp.

207 {
208  gpuErrchk(cudaMallocHost(h_data, Nbytes));
209 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31

◆ gpuCopyFromCPUToGPU()

void gpuCopyFromCPUToGPU ( void *  data,
void *  d_data,
size_t  Nbytes 
)

Definition at line 278 of file cuda_xmipp_utils.cpp.

279 {
280  gpuErrchk(cudaMemcpy(d_data, data, Nbytes, cudaMemcpyHostToDevice));
281 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31

◆ gpuCopyFromCPUToGPUStream()

void gpuCopyFromCPUToGPUStream ( void *  data,
void *  d_data,
size_t  Nbytes,
myStreamHandle myStream 
)

Definition at line 293 of file cuda_xmipp_utils.cpp.

294 {
295  cudaStream_t *stream = (cudaStream_t*) myStream.ptr;
296  gpuErrchk(cudaMemcpyAsync(d_data, data, Nbytes, cudaMemcpyHostToDevice, *stream));
297 
298  //gpuErrchk(cudaStreamSynchronize(*stream));
299 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31

◆ gpuCopyFromGPUToCPU()

void gpuCopyFromGPUToCPU ( void *  d_data,
void *  data,
size_t  Nbytes 
)

Definition at line 283 of file cuda_xmipp_utils.cpp.

284 {
285  gpuErrchk(cudaMemcpy(data, d_data, Nbytes, cudaMemcpyDeviceToHost));
286 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31

◆ gpuCopyFromGPUToCPUStream()

void gpuCopyFromGPUToCPUStream ( void *  d_data,
void *  data,
size_t  Nbytes,
myStreamHandle myStream 
)

Definition at line 301 of file cuda_xmipp_utils.cpp.

302 {
303  cudaStream_t *stream = (cudaStream_t*) myStream.ptr;
304  gpuErrchk(cudaMemcpyAsync(data, d_data, Nbytes, cudaMemcpyDeviceToHost, *stream));
305 
306  gpuErrchk(cudaStreamSynchronize(*stream));
307  //cudaDeviceSynchronize();
308 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31

◆ gpuCopyFromGPUToGPU()

void gpuCopyFromGPUToGPU ( void *  d_dataFrom,
void *  d_dataTo,
size_t  Nbytes 
)

Definition at line 288 of file cuda_xmipp_utils.cpp.

289 {
290  gpuErrchk(cudaMemcpy(d_dataTo, d_dataFrom, Nbytes, cudaMemcpyDeviceToDevice));
291 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31

◆ gpuCopyFromGPUToGPUStream()

void gpuCopyFromGPUToGPUStream ( void *  d_dataFrom,
void *  d_dataTo,
size_t  Nbytes,
myStreamHandle myStream 
)

Definition at line 310 of file cuda_xmipp_utils.cpp.

311 {
312  cudaStream_t *stream = (cudaStream_t*) myStream.ptr;
313  gpuErrchk(cudaMemcpyAsync(d_dataTo, d_dataFrom, Nbytes, cudaMemcpyDeviceToDevice, *stream));
314 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31

◆ gpuFree()

void gpuFree ( void *  d_data)

Definition at line 201 of file cuda_xmipp_utils.cpp.

202 {
203  gpuErrchk(cudaFree(d_data));
204 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31

◆ gpuMalloc()

void gpuMalloc ( void **  d_data,
size_t  Nbytes 
)

Definition at line 196 of file cuda_xmipp_utils.cpp.

197 {
198  gpuErrchk(cudaMalloc(d_data, Nbytes));
199 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31

◆ gridFromBlock()

int gridFromBlock ( int  tasks,
int  Nthreads 
)

Definition at line 316 of file cuda_xmipp_utils.cpp.

317 {
318  int numBlk = tasks/Nthreads;
319  if(tasks%Nthreads>0)
320  numBlk++;
321  return numBlk;
322 }

◆ initializeIdentity()

void initializeIdentity ( float *  d_data,
float *  h_data,
int  Ndim,
myStreamHandle myStream 
)

Definition at line 216 of file cuda_xmipp_utils.cpp.

217 {
218  //float *matrices = new float[Ndim*9];
219  for(int i=0; i<Ndim; i++){
220  float aux_matrix[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
221  int offset=i*9;
222  for (int j=0; j<9; j++)
223  h_data[offset+j] = aux_matrix[j];
224  }
225  cudaStream_t *stream = (cudaStream_t*) myStream.ptr;
226  gpuErrchk(cudaMemcpyAsync((void*)d_data, h_data, Ndim*9*sizeof(float), cudaMemcpyHostToDevice, *stream));
227  //delete []matrices;
228 
229 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31
#define i
#define j

◆ mycufftDestroy()

void mycufftDestroy ( void *  ptr)

Definition at line 72 of file cuda_xmipp_utils.cpp.

73 {
74  cufftHandle *planPtr = (cufftHandle *)ptr;
75  cufftDestroy(*planPtr);
76  delete planPtr;
77 }
int cufftHandle
Definition: cuda_fft.h:41

◆ myStreamCreate()

void myStreamCreate ( myStreamHandle myStream)

Definition at line 63 of file cuda_xmipp_utils.cpp.

64 {
65  cudaStream_t *streamPtr = new cudaStream_t;
66  gpuErrchk(cudaStreamCreate(streamPtr));
67  myStream.ptr = (void*)streamPtr;
68  //printf("ptr %p\n", myStream.ptr);
69  //printf("streamPtr %p\n", streamPtr);
70 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31

◆ myStreamDestroy()

void myStreamDestroy ( void *  ptr)

Definition at line 57 of file cuda_xmipp_utils.cpp.

58 {
59  cudaStream_t *streamPtr = (cudaStream_t *)ptr;
60  cudaStreamDestroy(*streamPtr);
61 }

◆ setRotationMatrix()

void setRotationMatrix ( float *  d_data,
float *  ang,
int  Ndim,
myStreamHandle myStream 
)

Definition at line 260 of file cuda_xmipp_utils.cpp.

261 {
262 
263  float *rad_vector;
264  gpuErrchk(cudaMallocHost((void**)&rad_vector, sizeof(float)*Ndim*9));
265 
266  for(int i=0; i<Ndim; i++){
267  float rad = (float)(-ang[i]*M_PI/180);
268  float matrix[9] = {cosf(rad), -sinf(rad), 0, sinf(rad), cosf(rad), 0, 0, 0, 1};
269  int offset=i*9;
270  for (int j=0; j<9; j++)
271  rad_vector[offset+j] = matrix[j];
272  }
273  cudaStream_t *stream = (cudaStream_t*) myStream.ptr;
274  gpuErrchk(cudaMemcpyAsync((void*)d_data, rad_vector, Ndim*9*sizeof(float), cudaMemcpyHostToDevice, *stream));
275  delete []rad_vector;
276 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31
#define i
#define j

◆ setTranslationMatrix()

void setTranslationMatrix ( float *  d_data,
float *  posX,
float *  posY,
int  Ndim,
myStreamHandle myStream 
)

Definition at line 237 of file cuda_xmipp_utils.cpp.

238 {
239  float *matrices;
240  gpuErrchk(cudaMallocHost((void**)&matrices, sizeof(float)*Ndim*9));
241 
242  for(int i=0; i<Ndim; i++){
243  float aux_matrix[9] = {1, 0, -posX[i], 0, 1, -posY[i], 0, 0, 1};
244  int offset=i*9;
245  //memcpy(&matrices[offset], &aux_matrix, 9*sizeof(float));
246  for (int j=0; j<9; j++)
247  matrices[offset+j] = aux_matrix[j];
248  }
249  cudaStream_t *stream = (cudaStream_t*) myStream.ptr;
250  gpuErrchk(cudaMemcpyAsync((void*)d_data, matrices, Ndim*9*sizeof(float), cudaMemcpyHostToDevice, *stream));
251  delete []matrices;
252 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31
#define i
#define j