Xmipp  v3.23.11-Nereus
Classes | Public Member Functions | Public Attributes | Friends | List of all members
VolumeHalvesRestorator< T > Class Template Reference

#include <cuda_volume_halves_restorator.h>

Classes

struct  Builder
 

Public Member Functions

 VolumeHalvesRestorator (int verbosity, unsigned denoisingIters, unsigned deconvolutionIters, T sigma, T lambda, unsigned differenceIters, T Kdiff, T bankStep, T bankOverlap, unsigned weightFun, T weightPower)
 
void apply (const MultidimArray< T > &volume1, const MultidimArray< T > &volume2, const int *mask)
 
const MultidimArray< T > & getReconstructedVolume1 () const
 
const MultidimArray< T > & getReconstructedVolume2 () const
 
const MultidimArray< T > & getFilterBankVolume () const
 
const MultidimArray< T > & getDeconvolvedS () const
 
const MultidimArray< T > & getConvolvedS () const
 
const MultidimArray< T > & getAverageDifference () const
 

Public Attributes

struct {
   T *   d_R2
 
   Complex *   d_fVol
 
   Complex *   d_fV1
 
   Complex *   d_fV2
 
   size_t   fourier_size
 
restorationPointers
 

Friends

std::ostream & operator<< (std::ostream &out, const VolumeHalvesRestorator &r)
 

Detailed Description

template<typename T>
class VolumeHalvesRestorator< T >

Definition at line 41 of file cuda_volume_halves_restorator.h.

Constructor & Destructor Documentation

◆ VolumeHalvesRestorator()

template<typename T>
VolumeHalvesRestorator< T >::VolumeHalvesRestorator ( int  verbosity,
unsigned  denoisingIters,
unsigned  deconvolutionIters,
sigma,
lambda,
unsigned  differenceIters,
Kdiff,
bankStep,
bankOverlap,
unsigned  weightFun,
weightPower 
)
inline

Definition at line 135 of file cuda_volume_halves_restorator.h.

137  : verbosity(verbosity)
138  , denoisingIters(denoisingIters)
139  , deconvolutionIters(deconvolutionIters)
140  , sigma(sigma)
141  , lambda(lambda)
142  , differenceIters(differenceIters)
143  , Kdiff(Kdiff)
144  , bankStep(bankStep)
145  , bankOverlap(bankOverlap)
146  , weightFun(weightFun)
147  , weightPower(weightPower) {}

Member Function Documentation

◆ apply()

template<typename T >
void VolumeHalvesRestorator< T >::apply ( const MultidimArray< T > &  volume1,
const MultidimArray< T > &  volume2,
const int *  mask 
)

Definition at line 345 of file cuda_volume_halves_restorator.cpp.

345  {
346  setSizes(volume1);
347  createFFTPlans();
348  initializeFilter();
349 
350  T* d_volume1;
351  T* d_volume2;
352  int* d_mask = nullptr;
353 
354  gpuErrchk( cudaMalloc(&d_volume1, memsize) );
355  gpuErrchk( cudaMalloc(&d_volume2, memsize) );
356  if (mask) {
357  gpuErrchk( cudaMalloc(&d_mask, volume_size * sizeof(int)) );
358  }
359 
360  gpuErrchk( cudaMemcpy(d_volume1, volume1.data, memsize, cudaMemcpyHostToDevice) );
361  gpuErrchk( cudaMemcpy(d_volume2, volume2.data, memsize, cudaMemcpyHostToDevice) );
362  if (mask) {
363  gpuErrchk( cudaMemcpy(d_mask, mask, volume_size * sizeof(int), cudaMemcpyHostToDevice) );
364  }
365 
366  gpuErrchk( cudaMalloc(&d_buf1, memsize) );
367  gpuErrchk( cudaMalloc(&d_buf2, memsize) );
368  gpuErrchk( cudaMalloc(&d_cbuf1, fourier_memsize) );
369  gpuErrchk( cudaMalloc(&d_cbuf2, fourier_memsize) );
370 
371  denoise(d_volume1, d_volume2, d_mask);
372  deconvolution(d_volume1, d_volume2);
373  filterBank(d_volume1, d_volume2);
374  difference(d_volume1, d_volume2, d_mask);
375 
376  reconstructedVolume1.resize(zdim, ydim, xdim);
377  gpuErrchk( cudaMemcpy(reconstructedVolume1.data, d_volume1, memsize, cudaMemcpyDeviceToHost) );
378  reconstructedVolume2.resize(zdim, ydim, xdim);
379  gpuErrchk( cudaMemcpy(reconstructedVolume2.data, d_volume2, memsize, cudaMemcpyDeviceToHost) );
380 
381  gpuErrchk( cudaFree(d_R2) );
382  gpuErrchk( cudaFree(d_mask) );
383  gpuErrchk( cudaFree(d_buf1) );
384  gpuErrchk( cudaFree(d_buf2) );
385  gpuErrchk( cudaFree(d_cbuf1) );
386  gpuErrchk( cudaFree(d_cbuf2) );
387  CudaFFT<T>::release(planForward);
388  CudaFFT<T>::release(planBackward);
389 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31
void release() final
Definition: cuda_fft.cpp:108

◆ getAverageDifference()

template<typename T>
const MultidimArray<T>& VolumeHalvesRestorator< T >::getAverageDifference ( ) const
inline

Definition at line 163 of file cuda_volume_halves_restorator.h.

163 { return averageDifference; }

◆ getConvolvedS()

template<typename T>
const MultidimArray<T>& VolumeHalvesRestorator< T >::getConvolvedS ( ) const
inline

Definition at line 162 of file cuda_volume_halves_restorator.h.

162 { return convolvedS; }

◆ getDeconvolvedS()

template<typename T>
const MultidimArray<T>& VolumeHalvesRestorator< T >::getDeconvolvedS ( ) const
inline

Definition at line 161 of file cuda_volume_halves_restorator.h.

161 { return deconvolvedS; }

◆ getFilterBankVolume()

template<typename T>
const MultidimArray<T>& VolumeHalvesRestorator< T >::getFilterBankVolume ( ) const
inline

Definition at line 160 of file cuda_volume_halves_restorator.h.

160 { return filterBankVolume; }

◆ getReconstructedVolume1()

template<typename T>
const MultidimArray<T>& VolumeHalvesRestorator< T >::getReconstructedVolume1 ( ) const
inline

Definition at line 158 of file cuda_volume_halves_restorator.h.

158 { return reconstructedVolume1; }

◆ getReconstructedVolume2()

template<typename T>
const MultidimArray<T>& VolumeHalvesRestorator< T >::getReconstructedVolume2 ( ) const
inline

Definition at line 159 of file cuda_volume_halves_restorator.h.

159 { return reconstructedVolume2; }

Friends And Related Function Documentation

◆ operator<<

template<typename T>
std::ostream& operator<< ( std::ostream &  out,
const VolumeHalvesRestorator< T > &  r 
)
friend

Definition at line 119 of file cuda_volume_halves_restorator.h.

119  {
120  out << "VolumeHalvesRestoration parameters:" << std::endl
121  << " Denoising Iterations:" << r.denoisingIters << std::endl
122  << " Deconvolution Iterations: " << r.deconvolutionIters << std::endl
123  << " Sigma0: " << r.sigma << std::endl
124  << " Lambda: " << r.lambda << std::endl
125  << " Bank step:" << r.bankStep << std::endl
126  << " Bank overlap:" << r.bankOverlap << std::endl
127  << " Weight fun:" << r.weightFun << std::endl
128  << " Weight power:" << r.weightPower << std::endl
129  << " Difference Iterations: " << r.differenceIters << std::endl
130  << " Kdiff: " << r.Kdiff << std::endl
131  ;
132  return out;
133  }

Member Data Documentation

◆ d_fV1

template<typename T>
Complex* VolumeHalvesRestorator< T >::d_fV1

Definition at line 172 of file cuda_volume_halves_restorator.h.

◆ d_fV2

template<typename T>
Complex* VolumeHalvesRestorator< T >::d_fV2

Definition at line 173 of file cuda_volume_halves_restorator.h.

◆ d_fVol

template<typename T>
Complex* VolumeHalvesRestorator< T >::d_fVol

Definition at line 171 of file cuda_volume_halves_restorator.h.

◆ restorationPointers

struct { ... } VolumeHalvesRestorator< T >::restorationPointers

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