Xmipp  v3.23.11-Nereus
Static Public Member Functions | List of all members
Gpu::VolumeRestorationKernels< T > Class Template Reference

#include <cuda_volume_restoration_kernels.h>

Static Public Member Functions

static void computeWeights (const T *d_Vfiltered1, const T *d_Vfiltered2, T *d_V1r, T *d_V2r, T *d_S, size_t volume_size, const Gpu::CDF< T > &cdf_mN, T weightPower, int weightFun)
 
static void filterFourierVolume (const T *d_R2, const std::complex< T > *d_fV, std::complex< T > *d_buffer, size_t volume_size, T w2, T w2Step)
 
static void computeAveragePositivity (const T *d_V1, const T *d_V2, T *d_S, size_t volume_size)
 
static void computeAveragePositivity (const T *d_V1, const T *d_V2, T *d_S, const int *d_mask, size_t volume_size)
 
static void filterS (const T *d_R2, std::complex< T > *d_fVol, size_t volume_size)
 
static void maskForCDF (T *__restrict__ d_aux, const T *__restrict__ d_S, const int *__restrict__ d_mask, size_t volume_size)
 
static void maskWithNoiseProbability (T *d_V, const Gpu::CDF< T > &cdf_S, const Gpu::CDF< T > &cdf_N, size_t volume_size)
 
static void deconvolveRestored (std::complex< T > *d_fVol, std::complex< T > *d_fV1, std::complex< T > *d_fV2, const T *d_R2, T K1, T K2, T lambda, size_t volume_size, size_t fourier_size)
 
static void convolveFourierVolume (std::complex< T > *d_fVol, const T *d_R2, T K, size_t volume_size)
 
static void normalizeForFFT (T *d_V1, T *d_V2, size_t volume_size)
 
static void normalizeForFFT (T *d_V1, size_t volume_size)
 
static void restorationSigmaCostError (T &error, const std::complex< T > *d_fVol, const std::complex< T > *d_fV1, const std::complex< T > *d_fV2, const T *__restrict__ d_R2, T K1, T K2, size_t fourier_size)
 
static void computeDiffAndAverage (const T *__restrict__ d_V1, const T *__restrict__ d_V2, T *__restrict__ d_S, T *__restrict__ d_N, size_t volume_size)
 
static std::pair< T, T > computeAvgStd (const T *__restrict__ d_N, size_t volume_size)
 
static std::pair< T, T > computeAvgStdWithMask (const T *__restrict__ d_N, const int *__restrict__ d_mask, size_t mask_size, size_t volume_size)
 
static void computeDifference (T *__restrict__ d_V1, T *__restrict__ d_V2, const T *__restrict__ d_S, const T *__restrict__ d_N, T k, size_t volume_size)
 
static size_t computeMaskSize (const int *__restrict__ d_mask, size_t volume_size)
 
static void multiplyByConstant (T *__restrict__ d_array, T c, size_t volume_size)
 

Detailed Description

template<typename T>
class Gpu::VolumeRestorationKernels< T >

Definition at line 43 of file cuda_volume_restoration_kernels.h.

Member Function Documentation

◆ computeAveragePositivity() [1/2]

template<typename T >
void Gpu::VolumeRestorationKernels< T >::computeAveragePositivity ( const T *  d_V1,
const T *  d_V2,
T *  d_S,
size_t  volume_size 
)
static

Definition at line 69 of file cuda_volume_restoration_kernels.cpp.

69  {
70  unsigned block_size = 256;
71 
72  dim3 dimBlock{ block_size };
73  dim3 dimGrid{ static_cast<unsigned>( ceil( static_cast<double>(volume_size) / dimBlock.x) ) };
74 
75  computeAveragePositivityKernel<<<dimGrid, dimBlock>>>(d_V1, d_V2, d_S, volume_size, static_cast<T>(1.0) / volume_size);
76  gpuErrchk(cudaPeekAtLastError());
77 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31

◆ computeAveragePositivity() [2/2]

template<typename T >
void Gpu::VolumeRestorationKernels< T >::computeAveragePositivity ( const T *  d_V1,
const T *  d_V2,
T *  d_S,
const int *  d_mask,
size_t  volume_size 
)
static

Definition at line 80 of file cuda_volume_restoration_kernels.cpp.

80  {
81  unsigned block_size = 256;
82 
83  dim3 dimBlock{ block_size };
84  dim3 dimGrid{ static_cast<unsigned>( ceil( static_cast<double>(volume_size) / dimBlock.x) ) };
85 
86  computeAveragePositivityKernel<<<dimGrid, dimBlock>>>(d_V1, d_V2, d_S, d_mask, volume_size, static_cast<T>(1.0) / volume_size);
87  gpuErrchk(cudaPeekAtLastError());
88 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31

◆ computeAvgStd()

template<typename T >
std::pair< T, T > Gpu::VolumeRestorationKernels< T >::computeAvgStd ( const T *__restrict__  d_N,
size_t  volume_size 
)
static

Definition at line 222 of file cuda_volume_restoration_kernels.cpp.

222  {
223  const T avg = thrust::reduce(thrust::device, d_N, d_N + volume_size);
224 
225  auto square_kernel = [=] __device__ (T x) {
226  return x * x;
227  };
228 
229  const T std = thrust::transform_reduce(thrust::device, d_N, d_N + volume_size, square_kernel, static_cast<T>(0), thrust::plus<T>());
230 
231  return normAvgStd(avg, std, volume_size);
232 }
doublereal * x
std::pair< T, T > normAvgStd(T avg, T std, size_t size)

◆ computeAvgStdWithMask()

template<typename T >
std::pair< T, T > Gpu::VolumeRestorationKernels< T >::computeAvgStdWithMask ( const T *__restrict__  d_N,
const int *__restrict__  d_mask,
size_t  mask_size,
size_t  volume_size 
)
static

Definition at line 235 of file cuda_volume_restoration_kernels.cpp.

235  {
236  auto masked_k = [=] __device__ (int n) {
237  if (d_mask[n]) {
238  return d_N[n];
239  } else {
240  return static_cast<T>(0);
241  }
242  };
243 
244  const T avg = thrust::transform_reduce(thrust::device, thrust::counting_iterator<int>(0), thrust::counting_iterator<int>(volume_size), masked_k,
245  static_cast<T>(0), thrust::plus<T>());
246 
247  auto masked_square_k = [=] __device__ (int n) {
248  if (d_mask[n]) {
249  return d_N[n] * d_N[n];
250  } else {
251  return static_cast<T>(0);
252  }
253  };
254 
255  const T std = thrust::transform_reduce(thrust::device, thrust::counting_iterator<int>(0), thrust::counting_iterator<int>(volume_size), masked_square_k,
256  static_cast<T>(0), thrust::plus<T>());
257 
258  return normAvgStd(avg, std, mask_size);
259 }
std::pair< T, T > normAvgStd(T avg, T std, size_t size)
int * n

◆ computeDiffAndAverage()

template<typename T >
void Gpu::VolumeRestorationKernels< T >::computeDiffAndAverage ( const T *__restrict__  d_V1,
const T *__restrict__  d_V2,
T *__restrict__  d_S,
T *__restrict__  d_N,
size_t  volume_size 
)
static

Definition at line 197 of file cuda_volume_restoration_kernels.cpp.

197  {
198  auto k = [=] __device__ (int n) {
199  d_N[n] = d_V1[n] - d_V2[n];
200  d_S[n] = (d_V1[n] + d_V2[n]) * static_cast<T>(0.5);
201  };
202 
203  thrust::for_each_n(thrust::device, thrust::counting_iterator<int>(0), volume_size, k);
204 }
ql0001_ & k(htemp+1),(cvec+1),(atemp+1),(bj+1),(bl+1),(bu+1),(x+1),(clamda+1), &iout, infoqp, &zero,(w+1), &lenw,(iw+1), &leniw, &glob_grd.epsmac
int * n

◆ computeDifference()

template<typename T >
void Gpu::VolumeRestorationKernels< T >::computeDifference ( T *__restrict__  d_V1,
T *__restrict__  d_V2,
const T *__restrict__  d_S,
const T *__restrict__  d_N,
k,
size_t  volume_size 
)
static

Definition at line 262 of file cuda_volume_restoration_kernels.cpp.

262  {
263  auto ker = [=] __device__ (int n) {
264  const T Nn = d_N[n];
265  const T w = exp(k * Nn * Nn);
266  const T s = d_S[n];
267  d_V1[n] = s + (d_V1[n] - s) * w;
268  d_V2[n] = s + (d_V2[n] - s) * w;
269  };
270 
271  thrust::for_each_n(thrust::device, thrust::counting_iterator<int>(0), volume_size, ker);
272 }
doublereal * w
ql0001_ & k(htemp+1),(cvec+1),(atemp+1),(bj+1),(bl+1),(bu+1),(x+1),(clamda+1), &iout, infoqp, &zero,(w+1), &lenw,(iw+1), &leniw, &glob_grd.epsmac
int * n

◆ computeMaskSize()

template<typename T >
size_t Gpu::VolumeRestorationKernels< T >::computeMaskSize ( const int *__restrict__  d_mask,
size_t  volume_size 
)
static

Definition at line 275 of file cuda_volume_restoration_kernels.cpp.

275  {
276  return thrust::reduce(thrust::device, d_mask, d_mask + volume_size);
277 }

◆ computeWeights()

template<typename T >
void Gpu::VolumeRestorationKernels< T >::computeWeights ( const T *  d_Vfiltered1,
const T *  d_Vfiltered2,
T *  d_V1r,
T *  d_V2r,
T *  d_S,
size_t  volume_size,
const Gpu::CDF< T > &  cdf_mN,
weightPower,
int  weightFun 
)
static

Definition at line 44 of file cuda_volume_restoration_kernels.cpp.

45  {
46 
47  unsigned block_size = 256;
48 
49  dim3 dimBlock{ block_size };
50  dim3 dimGrid{ static_cast<unsigned>( ceil( static_cast<double>(volume_size) / dimBlock.x) ) };
51 
52  computeWeightsKernel<<<dimGrid, dimBlock>>>(d_Vfiltered1, d_Vfiltered2, d_V1r, d_V2r, d_S, cdf_mN.d_x,
53  cdf_mN.d_probXLessThanx, cdf_mN.d_V, volume_size, cdf_mN.Nsteps, weightPower, weightFun);
54  gpuErrchk(cudaPeekAtLastError());
55 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31
T * d_x
Definition: cuda_cdf.h:37
T * d_probXLessThanx
Definition: cuda_cdf.h:38
T * d_V
Definition: cuda_cdf.h:36
T Nsteps
Definition: cuda_cdf.h:43

◆ convolveFourierVolume()

template<typename T >
void Gpu::VolumeRestorationKernels< T >::convolveFourierVolume ( std::complex< T > *  d_fVol,
const T *  d_R2,
K,
size_t  volume_size 
)
static

Definition at line 135 of file cuda_volume_restoration_kernels.cpp.

135  {
136  unsigned block_size = 256;
137 
138  dim3 dimBlock{ block_size };
139  dim3 dimGrid{ static_cast<unsigned>( ceil( static_cast<double>(volume_size) / dimBlock.x) ) };
140 
141  convolveFourierVolumeKernel<<<dimGrid, dimBlock>>>((vec2_type<T>*)d_fVol, d_R2, K, volume_size);
142  gpuErrchk(cudaPeekAtLastError());
143 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31
typename Vec2< T >::type vec2_type
Definition: cuda_vec2.h:20
constexpr int K

◆ deconvolveRestored()

template<typename T >
void Gpu::VolumeRestorationKernels< T >::deconvolveRestored ( std::complex< T > *  d_fVol,
std::complex< T > *  d_fV1,
std::complex< T > *  d_fV2,
const T *  d_R2,
K1,
K2,
lambda,
size_t  volume_size,
size_t  fourier_size 
)
static

Definition at line 124 of file cuda_volume_restoration_kernels.cpp.

124  {
125  unsigned block_size = 256;
126 
127  dim3 dimBlock{ block_size };
128  dim3 dimGrid{ static_cast<unsigned>( ceil( static_cast<double>(fourier_size) / dimBlock.x) ) };
129 
130  deconvolveRestoredKernel<<<dimGrid, dimBlock>>>((vec2_type<T>*)d_fVol, (vec2_type<T>*)d_fV1, (vec2_type<T>*)d_fV2, d_R2, K1, K2, lambda, fourier_size, static_cast<T>(1.0) / volume_size);
131  gpuErrchk(cudaPeekAtLastError());
132 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31
double * lambda
typename Vec2< T >::type vec2_type
Definition: cuda_vec2.h:20

◆ filterFourierVolume()

template<typename T >
void Gpu::VolumeRestorationKernels< T >::filterFourierVolume ( const T *  d_R2,
const std::complex< T > *  d_fV,
std::complex< T > *  d_buffer,
size_t  volume_size,
w2,
w2Step 
)
static

Definition at line 58 of file cuda_volume_restoration_kernels.cpp.

58  {
59  unsigned block_size = 256;
60 
61  dim3 dimBlock{ block_size };
62  dim3 dimGrid{ static_cast<unsigned>( ceil( static_cast<double>(volume_size) / dimBlock.x) ) };
63 
64  filterFourierVolumeKernel<<<dimGrid, dimBlock>>>(d_R2, (vec2_type<T>*)d_fV, (vec2_type<T>*)d_buffer, volume_size, w2, w2Step);
65  gpuErrchk(cudaPeekAtLastError());
66 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31
typename Vec2< T >::type vec2_type
Definition: cuda_vec2.h:20

◆ filterS()

template<typename T >
void Gpu::VolumeRestorationKernels< T >::filterS ( const T *  d_R2,
std::complex< T > *  d_fVol,
size_t  volume_size 
)
static

Definition at line 91 of file cuda_volume_restoration_kernels.cpp.

91  {
92  unsigned block_size = 256;
93 
94  dim3 dimBlock{ block_size };
95  dim3 dimGrid{ static_cast<unsigned>( ceil( static_cast<double>(volume_size) / dimBlock.x) ) };
96 
97  filterSKernel<<<dimGrid, dimBlock>>>(d_R2, (vec2_type<T>*)d_fVol, volume_size);
98  gpuErrchk(cudaPeekAtLastError());
99 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31
typename Vec2< T >::type vec2_type
Definition: cuda_vec2.h:20

◆ maskForCDF()

template<typename T >
void Gpu::VolumeRestorationKernels< T >::maskForCDF ( T *__restrict__  d_aux,
const T *__restrict__  d_S,
const int *__restrict__  d_mask,
size_t  volume_size 
)
static

Definition at line 102 of file cuda_volume_restoration_kernels.cpp.

102  {
103  auto k = [] __device__ (int x) {
104  return x;
105  };
106 
107  thrust::copy_if(thrust::device, d_S, d_S + volume_size, d_mask, d_aux, k);
108 }
doublereal * x
ql0001_ & k(htemp+1),(cvec+1),(atemp+1),(bj+1),(bl+1),(bu+1),(x+1),(clamda+1), &iout, infoqp, &zero,(w+1), &lenw,(iw+1), &leniw, &glob_grd.epsmac

◆ maskWithNoiseProbability()

template<typename T >
void Gpu::VolumeRestorationKernels< T >::maskWithNoiseProbability ( T *  d_V,
const Gpu::CDF< T > &  cdf_S,
const Gpu::CDF< T > &  cdf_N,
size_t  volume_size 
)
static

Definition at line 111 of file cuda_volume_restoration_kernels.cpp.

111  {
112  unsigned block_size = 256;
113 
114  dim3 dimBlock{ block_size };
115  dim3 dimGrid{ static_cast<unsigned>( ceil( static_cast<double>(volume_size) / dimBlock.x) ) };
116 
117  maskWithNoiseProbabilityKernel<<<dimGrid, dimBlock>>>(d_V, cdf_S.d_x,
118  cdf_S.d_probXLessThanx, cdf_S.d_V, cdf_S.Nsteps, cdf_S.volume_size, cdf_N.d_x,
119  cdf_N.d_probXLessThanx, cdf_N.d_V, cdf_N.Nsteps, cdf_N.volume_size);
120  gpuErrchk(cudaPeekAtLastError());
121 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31
T * d_x
Definition: cuda_cdf.h:37
size_t volume_size
Definition: cuda_cdf.h:40
T * d_probXLessThanx
Definition: cuda_cdf.h:38
T * d_V
Definition: cuda_cdf.h:36
T Nsteps
Definition: cuda_cdf.h:43

◆ multiplyByConstant()

template<typename T >
void Gpu::VolumeRestorationKernels< T >::multiplyByConstant ( T *__restrict__  d_array,
c,
size_t  volume_size 
)
static

Definition at line 280 of file cuda_volume_restoration_kernels.cpp.

280  {
281  auto k = [=] __device__ (int n) {
282  d_array[n] = d_array[n] * c;
283  };
284 
285  thrust::for_each_n(thrust::device, thrust::counting_iterator<int>(0), volume_size, k);
286 }
doublereal * c
ql0001_ & k(htemp+1),(cvec+1),(atemp+1),(bj+1),(bl+1),(bu+1),(x+1),(clamda+1), &iout, infoqp, &zero,(w+1), &lenw,(iw+1), &leniw, &glob_grd.epsmac
int * n

◆ normalizeForFFT() [1/2]

template<typename T >
void Gpu::VolumeRestorationKernels< T >::normalizeForFFT ( T *  d_V1,
T *  d_V2,
size_t  volume_size 
)
static

Definition at line 146 of file cuda_volume_restoration_kernels.cpp.

146  {
147  unsigned block_size = 256;
148 
149  dim3 dimBlock{ block_size };
150  dim3 dimGrid{ static_cast<unsigned>( ceil( static_cast<double>(volume_size) / dimBlock.x) ) };
151 
152  normalizeForFFTKernel<<<dimGrid, dimBlock>>>(d_V1, d_V2, volume_size, static_cast<T>(1.0) / volume_size);
153  gpuErrchk(cudaPeekAtLastError());
154 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31

◆ normalizeForFFT() [2/2]

template<typename T >
void Gpu::VolumeRestorationKernels< T >::normalizeForFFT ( T *  d_V1,
size_t  volume_size 
)
static

Definition at line 157 of file cuda_volume_restoration_kernels.cpp.

157  {
158  unsigned block_size = 256;
159 
160  dim3 dimBlock{ block_size };
161  dim3 dimGrid{ static_cast<unsigned>( ceil( static_cast<double>(volume_size) / dimBlock.x) ) };
162 
163  normalizeForFFTKernel<<<dimGrid, dimBlock>>>(d_V1, volume_size, static_cast<T>(1.0) / volume_size);
164  gpuErrchk(cudaPeekAtLastError());
165 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31

◆ restorationSigmaCostError()

template<typename T >
void Gpu::VolumeRestorationKernels< T >::restorationSigmaCostError ( T &  error,
const std::complex< T > *  d_fVol,
const std::complex< T > *  d_fV1,
const std::complex< T > *  d_fV2,
const T *__restrict__  d_R2,
K1,
K2,
size_t  fourier_size 
)
static

Definition at line 168 of file cuda_volume_restoration_kernels.cpp.

168  {
169  const T inv_size = 1.0 / (2 * fourier_size);
170 
171  const vec2_type<T>* __restrict__ d_fVol = (vec2_type<T>*)_d_fVol;
172  const vec2_type<T>* __restrict__ d_fV1 = (vec2_type<T>*)_d_fV1;
173  const vec2_type<T>* __restrict__ d_fV2 = (vec2_type<T>*)_d_fV2;
174 
175  auto error_func = [=] __device__ (int n) {
176  const T R2n = d_R2[n];
177  if (R2n <= 0.25) {
178  const T H1 = exp(K1 * R2n);
179  const T H2 = exp(K2 * R2n);
180 
181  const T diff1_x = (d_fVol[n].x*H1 - d_fV1[n].x) * inv_size;
182  const T diff1_y = (d_fVol[n].y*H1 - d_fV1[n].y) * inv_size;
183 
184  const T diff2_x = (d_fVol[n].x*H2 - d_fV2[n].x) * inv_size;
185  const T diff2_y = (d_fVol[n].y*H2 - d_fV2[n].y) * inv_size;
186 
187  return sqrt(diff1_x*diff1_x + diff1_y*diff1_y) + sqrt(diff2_x*diff2_x + diff2_y*diff2_y);
188  }
189 
190  return static_cast<T>(0.0);
191  };
192 
193  error = thrust::transform_reduce(thrust::device, thrust::counting_iterator<int>(0), thrust::counting_iterator<int>(fourier_size), error_func, static_cast<T>(0), thrust::plus<T>());
194 }
void sqrt(Image< double > &op)
typename Vec2< T >::type vec2_type
Definition: cuda_vec2.h:20
int * n

The documentation for this class was generated from the following files: