Xmipp  v3.23.11-Nereus
Public Member Functions | Protected Member Functions | List of all members
GeoTransformer< T > Class Template Reference

#include <cuda_gpu_geo_transformer.h>

Inheritance diagram for GeoTransformer< T >:
Inheritance graph
[legend]

Public Member Functions

 GeoTransformer ()
 
 ~GeoTransformer ()
 
void initForMatrix (size_t x, size_t y, size_t z)
 
void initLazyForMatrix (size_t x, size_t y=1, size_t z=1)
 
void initForBSpline (size_t inX, size_t inY, size_t inN, size_t splineX, size_t splineY, size_t splineN, const GPU &gpu)
 
void initLazyForBSpline (size_t inX, size_t inY, size_t inN, size_t splineX, size_t splineY, size_t splineN, const GPU &gpu)
 
void release ()
 
template<typename T_MAT >
void applyGeometry (int splineDegree, MultidimArray< T > &output, const MultidimArray< T > &input, const Matrix2D< T_MAT > &transform, bool isInv, bool wrap, T outside=0, const MultidimArray< T > *bCoeffsPtr=NULL)
 
void applyBSplineTransform (int splineDegree, MultidimArray< T > &output, const MultidimArray< T > &input, const std::pair< Matrix1D< T >, Matrix1D< T >> &coeffs, size_t imageIdx, T outside=0)
 
void test ()
 

Protected Member Functions

void applyBSplineTransformRef (int splineDegree, MultidimArray< T > &output, const MultidimArray< T > &input, const std::pair< Matrix1D< T >, Matrix1D< T >> &coeffs, size_t imageIdx, T outside=0)
 
void produceAndLoadCoeffs (const MultidimArray< T > &input)
 
std::unique_ptr< T[]> copy_out_d_in (size_t size) const
 

Detailed Description

template<typename T>
class GeoTransformer< T >

Definition at line 55 of file cuda_gpu_geo_transformer.h.

Constructor & Destructor Documentation

◆ GeoTransformer()

template<typename T>
GeoTransformer< T >::GeoTransformer ( )
inline

Constructor

Definition at line 59 of file cuda_gpu_geo_transformer.h.

59 { setDefaultValues(); };

◆ ~GeoTransformer()

template<typename T>
GeoTransformer< T >::~GeoTransformer ( )
inline

Definition at line 61 of file cuda_gpu_geo_transformer.h.

61  {
62  release();
63  }

Member Function Documentation

◆ applyBSplineTransform()

template<typename T >
void GeoTransformer< T >::applyBSplineTransform ( int  splineDegree,
MultidimArray< T > &  output,
const MultidimArray< T > &  input,
const std::pair< Matrix1D< T >, Matrix1D< T >> &  coeffs,
size_t  imageIdx,
outside = 0 
)

Apply local transformation defined by a BSpline

Parameters
splineDegreeused for interpolation
outputwhere resulting image will be stored
inputto process
coeffsfor the X and Y dimension of the input
imageIdxindex of the current image. This function assumes that multiple calls will be done and that interpolation is done also over time
outsidevalue of the output, where the interpolation does not store anything

Definition at line 187 of file cuda_gpu_geo_transformer.cpp.

190  {
191  checkRestrictions(3, output, input, coeffs, imageIdx);
192  auto stream = *(cudaStream_t*)gpu->stream();
193 
194  setOutputSize(output);
195  if ( splineDegree > 1 ) {
196  produceAndLoadCoeffs(input);
197  }
198  else
199  {
200  gpuErrchk(cudaMemcpyAsync(d_in, input.data, input.yxdim * sizeof(T), cudaMemcpyHostToDevice, stream));
201  }
202 
203  loadCoefficients(coeffs.first, coeffs.second);
204 
205  dim3 dimBlock(16, 16);
206  dim3 dimGrid(ceil(inX / (T) dimBlock.x), ceil((inY / (T) dimBlock.y) / (T) pixelsPerThread)); //more pixels
207 
208  // take into account end points
209  T hX = (splineX == 3) ? inX : (inX / (T) ((splineX - 3)));
210  T hY = (splineY == 3) ? inY : (inY / (T) ((splineY - 3)));
211  T hT = (splineN == 3) ? inN : (inN / (T) ((splineN - 3)));
212  T tPos = imageIdx / hT;
213 
214  switch (splineDegree) {
215  case 1:
216  applyLocalShiftGeometryKernelMorePixels<T, 1, pixelsPerThread><<<dimGrid, dimBlock, 0, stream>>>(d_coeffsX, d_coeffsY,
217  d_out, (int)inX, (int)inY, (int)inN,
218  d_in, imageIdx, (int)splineX, (int)splineY, (int)splineN,
219  hX, hY, tPos);
220  gpuErrchk(cudaPeekAtLastError());
221  break;
222  case 3:
223  applyLocalShiftGeometryKernelMorePixels<T, 3, pixelsPerThread><<<dimGrid, dimBlock, 0, stream>>>(d_coeffsX, d_coeffsY,
224  d_out, (int)inX, (int)inY, (int)inN,
225  d_in, imageIdx, (int)splineX, (int)splineY, (int)splineN,
226  hX, hY, tPos);
227  gpuErrchk(cudaPeekAtLastError());
228  break;
229  default:
230  REPORT_ERROR(ERR_NOT_IMPLEMENTED, formatString("applyBSplineTransform not implemented for spline degree %d.", splineDegree));
231  }
232 
233  gpuErrchk(
234  cudaMemcpyAsync(output.data, d_out, output.zyxdim * sizeof(T),
235  cudaMemcpyDeviceToHost, stream));
236 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31
void * stream() const
Definition: gpu.h:50
Case or algorithm not implemented yet.
Definition: xmipp_error.h:177
#define REPORT_ERROR(nerr, ErrormMsg)
Definition: xmipp_error.h:211
void produceAndLoadCoeffs(const MultidimArray< T > &input)
String formatString(const char *format,...)

◆ applyBSplineTransformRef()

template<typename T >
void GeoTransformer< T >::applyBSplineTransformRef ( int  splineDegree,
MultidimArray< T > &  output,
const MultidimArray< T > &  input,
const std::pair< Matrix1D< T >, Matrix1D< T >> &  coeffs,
size_t  imageIdx,
outside = 0 
)
protected

Definition at line 155 of file cuda_gpu_geo_transformer.cpp.

158  {
159  checkRestrictions(3, output, input, coeffs, imageIdx);
160 
161  loadOutput(output, outside);
162  produceAndLoadCoeffs(input);
163 
164  loadCoefficients(coeffs.first, coeffs.second);
165 
166  dim3 dimBlock(BLOCK_DIM_X, BLOCK_DIM_X);
167  dim3 dimGrid(ceil(inX / (T) dimBlock.x), ceil(inY / (T) dimBlock.y));
168 
169  switch (splineDegree) {
170  case 3:
171  applyLocalShiftGeometryKernel<T, 3><<<dimGrid, dimBlock>>>(d_coeffsX, d_coeffsY,
172  d_out, (int)inX, (int)inY, (int)inN,
173  d_in, imageIdx, (int)splineX, (int)splineY, (int)splineN);
174  gpuErrchk(cudaPeekAtLastError());
175  break;
176  default:
177  throw std::logic_error("not implemented");
178  }
179 
180  gpuErrchk(
181  cudaMemcpy(output.data, d_out, output.zyxdim * sizeof(T),
182  cudaMemcpyDeviceToHost));
183 
184 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31
void produceAndLoadCoeffs(const MultidimArray< T > &input)
#define BLOCK_DIM_X

◆ applyGeometry()

template<typename T >
template<typename T_MAT >
void GeoTransformer< T >::applyGeometry ( int  splineDegree,
MultidimArray< T > &  output,
const MultidimArray< T > &  input,
const Matrix2D< T_MAT > &  transform,
bool  isInv,
bool  wrap,
outside = 0,
const MultidimArray< T > *  bCoeffsPtr = NULL 
)

Apply geometry transformation Currently supported transformations: FIXME

  • 2D shift using cubic interpolation and mirrorOffBound
    Parameters
    splineDegreeused for interpolation
    outputwhere resulting image will be stored
    inputto process
    transformto apply
    isInvif true, 'transform' is expected to be from the resulting image to source image
    wraptrue to wrap after boundaris
    outsidevalue to be used when reading outside of the image and wrap == false
    bCoeffsPtrspline coefficients to use

Definition at line 240 of file cuda_gpu_geo_transformer.cpp.

243  {
244  checkRestrictions(splineDegree, output, input, transform);
245  if (transform.isIdentity()) {
246  typeCast(input, output);
247  return;
248  }
249 
250  loadTransform(transform, isInv);
251  loadOutput(output, outside);
252 
253  if (splineDegree > 1) {
254  if (NULL != bCoeffsPtr) {
255  loadInput(*bCoeffsPtr);
256  } else {
257  produceAndLoadCoeffs(input);
258  }
259  } else {
260  loadInput(input);
261  }
262 
263  if (input.getDim() == 2) {
264  if (wrap) {
265  applyGeometry_2D_wrap(splineDegree);
266  } else {
267  throw std::logic_error("Not implemented yet");
268  }
269  } else {
270  throw std::logic_error("Not implemented yet");
271  }
272 
273  gpuErrchk(
274  cudaMemcpy(output.data, d_out, output.zyxdim * sizeof(T),
275  cudaMemcpyDeviceToHost));
276 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31
bool isIdentity() const
Definition: matrix2d.cpp:1323
void produceAndLoadCoeffs(const MultidimArray< T > &input)
void typeCast(const Matrix1D< T1 > &v1, Matrix1D< T2 > &v2)
Definition: matrix1d.h:1227

◆ copy_out_d_in()

template<typename T >
std::unique_ptr< T[]> GeoTransformer< T >::copy_out_d_in ( size_t  size) const
protected

Definition at line 429 of file cuda_gpu_geo_transformer.cpp.

429  {
430  auto copy_in = std::unique_ptr<T[]>(new T[size]);
431 
432  cudaMemcpy(copy_in.get(), d_in , sizeof(T) * size, cudaMemcpyDeviceToHost);
433 
434  return copy_in;
435  }

◆ initForBSpline()

template<typename T >
void GeoTransformer< T >::initForBSpline ( size_t  inX,
size_t  inY,
size_t  inN,
size_t  splineX,
size_t  splineY,
size_t  splineN,
const GPU gpu 
)

Release previously obtained resources and initialize the transformer for processing images using BSpline coefficients. It also allocates all resources on GPU.

Parameters
sizesof the input images and number of images to be processed
numberof BSpline control points, including end points

Definition at line 71 of file cuda_gpu_geo_transformer.cpp.

73 {
74  release();
75 
76  this->inX = inX;
77  this->inY = inY;
78  this->inZ = 1;
79  this->inN = inN;
80  this->splineX = splineX;
81  this->splineY = splineY;
82  this->splineN = splineN;
83  // take into account end control points
84 
85  // padding for produceAndLoadCoeffs; Y dimension has to be a multiple of BLOCK_SIZE
86  const int BLOCK_SIZE = iirConvolve2D_Cardinal_BSpline_3_MirrorOffBoundKernels::BLOCK_SIZE;
87  const int Y_padded = (inY / BLOCK_SIZE) * BLOCK_SIZE + BLOCK_SIZE * (inY % BLOCK_SIZE != 0);
88 
89  size_t inOutSize = inX * inY;
90  size_t inOutSize_padded = inX* Y_padded;
91  size_t coeffsSize = splineX * splineY * splineN;
92  gpuErrchk(cudaMalloc((void** ) &d_coeffsX, coeffsSize * sizeof(T)));
93  gpuErrchk(cudaMalloc((void** ) &d_coeffsY, coeffsSize * sizeof(T)));
94  gpuErrchk(cudaMalloc((void** ) &d_in, inOutSize_padded * sizeof(T)));
95  gpuErrchk(cudaMalloc((void** ) &d_out, inOutSize * sizeof(T)));
96 
97  this->gpu = &gpu;
98 
99  isReadyForBspline = true;
100 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31

◆ initForMatrix()

template<typename T >
void GeoTransformer< T >::initForMatrix ( size_t  x,
size_t  y,
size_t  z 
)

Release previously obtained resources and initialize the transformer for processing images of given size. It also allocates all resources on GPU.

Parameters
xdim (inner-most) of the resulting image
ydim of the resulting image
zdim (outer-most) of the resulting image

Definition at line 49 of file cuda_gpu_geo_transformer.cpp.

49  {
50  release();
51 
52  inX = x;
53  inY = y;
54  inZ = z;
55  size_t matSize = (0 == z) ? 9 : 16;
56  gpuErrchk(cudaMalloc((void** ) &d_trInv, matSize * sizeof(T)));
57  gpuErrchk(cudaMalloc((void** ) &d_in, x * y * z * sizeof(T)));
58  gpuErrchk(cudaMalloc((void** ) &d_out, x * y * z * sizeof(T)));
59 
60  isReadyForMatrix = true;
61 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31
static double * y
doublereal * x
double z

◆ initLazyForBSpline()

template<typename T >
void GeoTransformer< T >::initLazyForBSpline ( size_t  inX,
size_t  inY,
size_t  inN,
size_t  splineX,
size_t  splineY,
size_t  splineN,
const GPU gpu 
)

Similar as the other init() function, except this method has no effect should the instance be already initialized. It is useful for example in a for loop, where first call will initialize resources and following calls will be ignored

Definition at line 103 of file cuda_gpu_geo_transformer.cpp.

105 {
106  if (!isReadyForBspline)
107  {
108  initForBSpline(inX, inY, inZ, splineX, splineY, splineN, gpu);
109  }
110 }
void initForBSpline(size_t inX, size_t inY, size_t inN, size_t splineX, size_t splineY, size_t splineN, const GPU &gpu)

◆ initLazyForMatrix()

template<typename T >
void GeoTransformer< T >::initLazyForMatrix ( size_t  x,
size_t  y = 1,
size_t  z = 1 
)

Similar as other init() function, except this method has no effect should the instance be already initialized. It is useful for example in a for loop, where first call will initialize resources and following calls will be ignored

Definition at line 64 of file cuda_gpu_geo_transformer.cpp.

64  {
65  if (!isReadyForMatrix) {
66  initForMatrix(x, y, z);
67  }
68 }
static double * y
void initForMatrix(size_t x, size_t y, size_t z)
doublereal * x
double z

◆ produceAndLoadCoeffs()

template<typename T >
void GeoTransformer< T >::produceAndLoadCoeffs ( const MultidimArray< T > &  input)
protected

Computes spline coefficients of the image and load them to GPU

Parameters
splineDegreeto be used
inputimage used to generate the coefficients

Definition at line 303 of file cuda_gpu_geo_transformer.cpp.

305 {
306  gpuErrchk(
307  cudaMemcpyAsync(d_in, input.data, input.yxdim * sizeof(T), cudaMemcpyHostToDevice, *(cudaStream_t*)gpu->stream()));
308 
309  iirConvolve2D_Cardinal_Bspline_3_MirrorOffBoundInplace(d_in, input.xdim, input.ydim, *(cudaStream_t*)gpu->stream());
310 }
#define gpuErrchk(code)
Definition: cuda_asserts.h:31
void * stream() const
Definition: gpu.h:50

◆ release()

template<typename T >
void GeoTransformer< T >::release ( )

Release all resources hold by this instance

Definition at line 32 of file cuda_gpu_geo_transformer.cpp.

32  {
33  cudaFree(d_in);
34  cudaFree(d_out);
35  cudaFree(d_trInv);
36  cudaFree(d_coeffsX);
37  cudaFree(d_coeffsY);
38  setDefaultValues();
39 }

◆ test()

template<typename T >
void GeoTransformer< T >::test ( )

Definition at line 113 of file cuda_gpu_geo_transformer.cpp.

113  {
114  Matrix1D<T> shift(2);
115  shift.vdata[0] = 0.45;
116  shift.vdata[1] = 0.62;
117  Matrix2D<T> transform;
118  translation2DMatrix(shift, transform, true);
119  test(transform);
120 }
template void translation2DMatrix(const Matrix1D< float > &, Matrix2D< float > &, bool inverse)
Definition: mask.h:36
Definition: ctf.h:38

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