Xmipp  v3.23.11-Nereus
cuda_gpu_correlation.h
Go to the documentation of this file.
1 
2 #ifndef CUDA_GPU_CORRELATION_H
3 #define CUDA_GPU_CORRELATION_H
4 #include <complex>
5 #include "cuda_xmipp_utils.h"
6 
8 void waitGpu(myStreamHandle &myStream, bool allStreams);
9 void calculateAbs (std::complex<float> *data, float *out, int size, myStreamHandle &myStream);
10 
15 public:
29 
30 };
31 
33 public:
42  int maskCount;
43  size_t Xdim, Ydim, XdimOrig, YdimOrig, XdimPolar, YdimPolar;
44  //AJ Xdim and Ydim are the padded sizes
45 
48 
50 
51  void produceSideInfo(mycufftHandle &myhandlePadded, mycufftHandle &myhandleMask, StructuresAux &myStructureAux, myStreamHandle &myStream);
52  void produceSideInfo(mycufftHandle &myhandlePadded, mycufftHandle &myhandleMask, StructuresAux &myStructureAux, GpuMultidimArrayAtGpu<float> &maskAutocorr, myStreamHandle &myStream);
53  void produceSideInfo(mycufftHandle &myhandlePadded, mycufftHandle &myhandleMask, StructuresAux &myStructureAux, GpuMultidimArrayAtGpu<float> &maskAutocorr, myStreamHandle &myStream, bool skip, mycufftHandle &ifftcb);
54 };
55 
56 
57 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);
58 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);
59 void apply_transform(GpuMultidimArrayAtGpu<float> &d_original_image, GpuMultidimArrayAtGpu<float> &d_transform_image, TransformMatrix<float> &transMat, myStreamHandle &myStream);
61 /*void calculateFFTPlanSize(mycufftHandle &myhandle);
62 void produceSideInfoCuda(GpuCorrelationAux &aux, mycufftHandle &myhandlePadded, mycufftHandle &myhandleMask, mycufftHandle &myhandlePolar, size_t memoryUsed);
63 void cuda_align_experimental_image(float *expImages, GpuCorrelationAux &d_referenceAux, GpuCorrelationAux &d_experimentalAux,
64  TransformMatrix<float> &transMat_tr, TransformMatrix<float> &transMat_rt, float *max_vector_tr, float *max_vector_rt,
65  int numImagesRef, GpuMultidimArrayAtGpu<float> &mask, double maxShift,
66  mycufftHandle &myhandlePadded, mycufftHandle &myhandleMask, mycufftHandle &myhandlePolar, size_t mdExpSize);
67 */
68 void cuda_calculate_correlation_two(GpuCorrelationAux &referenceAux, GpuCorrelationAux &experimentalAuxTR,
69  TransformMatrix<float> &transMatTR, float *max_vectorTR, int maxShift,
70  mycufftHandle &myhandlePaddedTR, bool mirror, StructuresAux &myStructureAuxTR,
71  myStreamHandle &myStreamTR,
72  GpuCorrelationAux &experimentalAuxRT, TransformMatrix<float> &transMatRT,
73  float *max_vectorRT, mycufftHandle &myhandlePaddedRT,
74  StructuresAux &myStructureAuxRT, myStreamHandle &myStreamRT,
76  mycufftHandle &ifftcb, bool saveMaxVector);
78 #endif
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)
GpuMultidimArrayAtGpu< float > polar2_gpu
void cuda_cart2polar(GpuMultidimArrayAtGpu< float > &image, GpuMultidimArrayAtGpu< float > &polar_image, GpuMultidimArrayAtGpu< float > &polar2_image, bool rotate, 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)
GpuMultidimArrayAtGpu< float > RefExpRealSpacePolar
void apply_transform(GpuMultidimArrayAtGpu< float > &d_original_image, GpuMultidimArrayAtGpu< float > &d_transform_image, TransformMatrix< float > &transMat, myStreamHandle &myStream)
GpuMultidimArrayAtGpu< std::complex< float > > d_projPolarSquaredFFT
void calculateAbs(std::complex< float > *data, float *out, int size, myStreamHandle &myStream)
GpuMultidimArrayAtGpu< float > d_pos_polar_max
void waitGpu(myStreamHandle &myStream, bool allStreams)
GpuMultidimArrayAtGpu< float > auxMax
GpuMultidimArrayAtGpu< float > RefExpRealSpace
GpuMultidimArrayAtGpu< float > padded_mask_gpu
GpuMultidimArrayAtGpu< float > padded_image2_gpu
#define rotate(a, i, j, k, l)
GpuMultidimArrayAtGpu< std::complex< float > > RefExpFourier
GpuMultidimArrayAtGpu< float > polar_gpu
GpuMultidimArrayAtGpu< float > maskAutocorrelation
GpuMultidimArrayAtGpu< std::complex< float > > d_projFFT
GpuMultidimArrayAtGpu< float > d_transform_image
GpuMultidimArrayAtGpu< std::complex< float > > RefExpFourierPolar
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)
GpuMultidimArrayAtGpu< float > d_out_max
GpuMultidimArrayAtGpu< std::complex< float > > d_projPolarFFT
GpuMultidimArrayAtGpu< float > d_original_image
GpuMultidimArrayAtGpu< float > d_NCCPolar
GpuMultidimArrayAtGpu< float > d_NCCPolar1D
GpuMultidimArrayAtGpu< std::complex< float > > MF2
GpuMultidimArrayAtGpu< float > d_denom
GpuMultidimArrayAtGpu< std::complex< float > > d_projSquaredFFT
GpuMultidimArrayAtGpu< float > d_pos_max
GpuMultidimArrayAtGpu< float > MFrealSpace
GpuMultidimArrayAtGpu< float > maxGpu
GpuMultidimArrayAtGpu< float > d_NCC
GpuMultidimArrayAtGpu< float > auxZero
GpuMultidimArrayAtGpu< float > d_out_polar_max
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)
GpuMultidimArrayAtGpu< std::complex< float > > MF
GpuMultidimArrayAtGpu< std::complex< float > > d_maskFFT
GpuMultidimArrayAtGpu< float > d_mask
GpuMultidimArrayAtGpu< float > padded_image_gpu