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

Go to the source code of this file.

Classes

class  StructuresAux
 
class  GpuCorrelationAux
 

Functions

void cuda_cart2polar (GpuMultidimArrayAtGpu< float > &image, GpuMultidimArrayAtGpu< float > &polar_image, GpuMultidimArrayAtGpu< float > &polar2_image, bool rotate, myStreamHandle &myStream)
 
void waitGpu (myStreamHandle &myStream, bool allStreams)
 
void calculateAbs (std::complex< float > *data, float *out, int size, myStreamHandle &myStream)
 
void cuda_calculate_correlation (GpuCorrelationAux &referenceAux, GpuCorrelationAux &experimentalAux, TransformMatrix< float > &transMat, float *max_vector, int maxShift, mycufftHandle &myhandlePadded, bool mirror, StructuresAux &myStructureAux, myStreamHandle &myStream, TransformMatrix< float > &resultTR, bool saveMaxVector)
 
void cuda_calculate_correlation_rotation (GpuCorrelationAux &referenceAux, GpuCorrelationAux &experimentalAux, TransformMatrix< float > &transMat, float *max_vector, int maxShift, mycufftHandle &myhandlePadded, bool mirror, StructuresAux &myStructureAux, myStreamHandle &myStream, TransformMatrix< float > &resultRT)
 
void apply_transform (GpuMultidimArrayAtGpu< float > &d_original_image, GpuMultidimArrayAtGpu< float > &d_transform_image, TransformMatrix< float > &transMat, myStreamHandle &myStream)
 
void padding_masking (GpuMultidimArrayAtGpu< float > &d_orig_image, GpuMultidimArrayAtGpu< float > &mask, GpuMultidimArrayAtGpu< float > &padded_image_gpu, GpuMultidimArrayAtGpu< float > &padded_image2_gpu, GpuMultidimArrayAtGpu< float > &padded_mask_gpu, bool experimental, myStreamHandle &myStream)
 
void cuda_calculate_correlation_two (GpuCorrelationAux &referenceAux, GpuCorrelationAux &experimentalAuxTR, TransformMatrix< float > &transMatTR, float *max_vectorTR, int maxShift, mycufftHandle &myhandlePaddedTR, bool mirror, StructuresAux &myStructureAuxTR, myStreamHandle &myStreamTR, GpuCorrelationAux &experimentalAuxRT, TransformMatrix< float > &transMatRT, float *max_vectorRT, mycufftHandle &myhandlePaddedRT, StructuresAux &myStructureAuxRT, myStreamHandle &myStreamRT, TransformMatrix< float > &resultTR, TransformMatrix< float > &resultRT, mycufftHandle &ifftcb, bool saveMaxVector)
 

Function Documentation

◆ calculateAbs()

void calculateAbs ( std::complex< float > *  data,
float *  out,
int  size,
myStreamHandle myStream 
)

Definition at line 1505 of file cuda_gpu_correlation.cpp.

1505  {
1506 
1507  cudaStream_t *stream = (cudaStream_t*) myStream.ptr;
1508  int numTh = 1024;
1509  int numBlk = size/numTh;
1510  if(size%numTh > 0)
1511  numBlk++;
1512  calcAbsKernel <<< numBlk, numTh, 0, *stream>>> ((cufftComplex*)data, out, size);
1513 
1514 }

◆ cuda_cart2polar()

void cuda_cart2polar ( GpuMultidimArrayAtGpu< float > &  image,
GpuMultidimArrayAtGpu< float > &  polar_image,
GpuMultidimArrayAtGpu< float > &  polar2_image,
bool  rotate,
myStreamHandle myStream 
)

Definition at line 1478 of file cuda_gpu_correlation.cpp.

1480 {
1481  cudaStream_t *stream = (cudaStream_t*) myStream.ptr;
1482  int numTh = 32;
1483  int numBlkx = polar_image.Xdim/numTh;
1484  if(polar_image.Xdim%numTh > 0)
1485  numBlkx++;
1486  int numBlky = polar_image.Ydim/numTh;
1487  if(polar_image.Ydim%numTh > 0)
1488  numBlky++;
1489 
1490  dim3 blockSize(numTh, numTh, 1);
1491  dim3 gridSize(numBlkx, numBlky, polar_image.Ndim);
1492 
1493  cart2polar <<< gridSize, blockSize, 0, *stream>>>
1494  (image.d_data, polar_image.d_data, polar2_image.d_data, polar_image.Ydim, polar_image.Xdim, polar_image.Ndim, image.Ydim, image.Xdim, rotate);
1495 }
#define rotate(a, i, j, k, l)

◆ waitGpu()

void waitGpu ( myStreamHandle myStream,
bool  allStreams 
)

Definition at line 1497 of file cuda_gpu_correlation.cpp.

1497  {
1498  if(!allStreams){
1499  cudaStream_t *stream = (cudaStream_t*) myStream.ptr;
1500  gpuErrchk(cudaStreamSynchronize(*stream));
1501  }else
1502  gpuErrchk(cudaDeviceSynchronize());
1503 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31