Grangeat-based 2D/3D image registration
Loading...
Searching...
No Matches
Functions
PyTorch Functions

Functions that have Python bindings generated and are accessible in the PyTorch extension. More...

Functions

at::Tensor reg23::GridSample3D_CPU (const at::Tensor &input, const at::Tensor &grid, const std::string &addressModeX, const std::string &addressModeY, const std::string &addressModeZ, c10::optional< at::Tensor > out)
 Sample the given 3D input tensor at the positions given in grid according to the given address mode using bilinear interpolation. This implementation is single-threaded.
 
__host__ at::Tensor reg23::GridSample3D_CUDA (const at::Tensor &input, const at::Tensor &grid, const std::string &addressModeX, const std::string &addressModeY, const std::string &addressModeZ, c10::optional< at::Tensor > out)
 An implementation of reg23::GridSample3D_CPU that uses CUDA parallelisation.
 
at::Tensor reg23::ProjectDRR_CPU (const at::Tensor &volume, const at::Tensor &voxelSpacing, const at::Tensor &homographyMatrixInverse, double sourceDistance, int64_t outputWidth, int64_t outputHeight, const at::Tensor &outputOffset, const at::Tensor &detectorSpacing)
 Generate a DRR from the given volume at the given transformation.
 
__host__ at::Tensor reg23::ProjectDRR_CUDA (const at::Tensor &volume, const at::Tensor &voxelSpacing, const at::Tensor &homographyMatrixInverse, double sourceDistance, int64_t outputWidth, int64_t outputHeight, const at::Tensor &outputOffset, const at::Tensor &detectorSpacing)
 An implementation of reg23::ProjectDRR_CPU that uses CUDA parallelisation.
 
at::Tensor reg23::ProjectDRR_backward_CPU (const at::Tensor &volume, const at::Tensor &voxelSpacing, const at::Tensor &homographyMatrixInverse, double sourceDistance, int64_t outputWidth, int64_t outputHeight, const at::Tensor &outputOffset, const at::Tensor &detectorSpacing, const at::Tensor &dLossDDRR)
 Evaluate the derivative of some scalar loss that is a function of a DRR projected from the given volume at the given transformation, with respect to the inverse homography matrix.
 
__host__ at::Tensor reg23::ProjectDRR_backward_CUDA (const at::Tensor &volume, const at::Tensor &voxelSpacing, const at::Tensor &homographyMatrixInverse, double sourceDistance, int64_t outputWidth, int64_t outputHeight, const at::Tensor &outputOffset, const at::Tensor &detectorSpacing, const at::Tensor &dLossDDRR)
 An implementation of reg23::ProjectDRR_backward_CPU that uses CUDA parallelisation.
 
__host__ at::Tensor reg23::ProjectDRRsBatched_CUDA (const at::Tensor &volume, const at::Tensor &voxelSpacing, const at::Tensor &invHMatrices, double sourceDistance, int64_t outputWidth, int64_t outputHeight, const at::Tensor &outputOffset, const at::Tensor &detectorSpacing)
 An implementation similar to reg23::ProjectDRR_CUDA that evaluates projections for multiple transformations in parallel.
 
at::Tensor reg23::ProjectDRRCuboidMask_CPU (const at::Tensor &volumeSize, const at::Tensor &voxelSpacing, const at::Tensor &homographyMatrixInverse, double sourceDistance, int64_t outputWidth, int64_t outputHeight, const at::Tensor &outputOffset, const at::Tensor &detectorSpacing)
 Calculate the distances of the intersections between the DRR-generation rays and the domain of the volume data.
 
__host__ at::Tensor reg23::ProjectDRRCuboidMask_CUDA (const at::Tensor &volumeSize, const at::Tensor &voxelSpacing, const at::Tensor &homographyMatrixInverse, double sourceDistance, int64_t outputWidth, int64_t outputHeight, const at::Tensor &outputOffset, const at::Tensor &detectorSpacing)
 An implementation of reg23::ProjectDRRCuboidMask_CPU that uses CUDA parallelisation.
 
at::Tensor reg23::Radon2D_CPU (const at::Tensor &image, const at::Tensor &imageSpacing, const at::Tensor &phiValues, const at::Tensor &rValues, int64_t samplesPerLine)
 Compute an approximation of the Radon transform of the given 2D image.
 
at::Tensor reg23::DRadon2DDR_CPU (const at::Tensor &image, const at::Tensor &imageSpacing, const at::Tensor &phiValues, const at::Tensor &rValues, int64_t samplesPerLine)
 Compute the derivative with respect to plane-origin distance of an approximation of the Radon transform of the given 2D image.
 
__host__ at::Tensor reg23::Radon2D_CUDA (const at::Tensor &image, const at::Tensor &imageSpacing, const at::Tensor &phiValues, const at::Tensor &rValues, int64_t samplesPerLine)
 An implementation of reg23::Radon2D_CPU that uses CUDA parallelisation.
 
__host__ at::Tensor reg23::DRadon2DDR_CUDA (const at::Tensor &image, const at::Tensor &imageSpacing, const at::Tensor &phiValues, const at::Tensor &rValues, int64_t samplesPerLine)
 An implementation of reg23::DRadon2DDR_CPU that uses CUDA parallelisation.
 
at::Tensor reg23::Radon3D_CPU (const at::Tensor &volume, const at::Tensor &volumeSpacing, const at::Tensor &phiValues, const at::Tensor &thetaValues, const at::Tensor &rValues, int64_t samplesPerDirection)
 Compute an approximation of the Radon transform of the given 3D volume.
 
at::Tensor reg23::DRadon3DDR_CPU (const at::Tensor &volume, const at::Tensor &volumeSpacing, const at::Tensor &phiValues, const at::Tensor &thetaValues, const at::Tensor &rValues, int64_t samplesPerDirection)
 Compute the derivative with respect to plane-origin distance of an approximation of the Radon transform of the given 3D volume.
 
__host__ at::Tensor reg23::Radon3D_CUDA (const at::Tensor &volume, const at::Tensor &volumeSpacing, const at::Tensor &phiValues, const at::Tensor &thetaValues, const at::Tensor &rValues, int64_t samplesPerDirection)
 An implementation of reg23::Radon3D_CPU that uses CUDA parallelisation.
 
__host__ at::Tensor reg23::DRadon3DDR_CUDA (const at::Tensor &volume, const at::Tensor &volumeSpacing, const at::Tensor &phiValues, const at::Tensor &thetaValues, const at::Tensor &rValues, int64_t samplesPerDirection)
 An implementation of reg23::DRadon3DDR_CPU that uses CUDA parallelisation.
 
at::Tensor reg23::ResampleSinogram3D_CPU (const at::Tensor &sinogram3d, const std::string &sinogramType, double rSpacing, const at::Tensor &projectionMatrix, const at::Tensor &phiValues, const at::Tensor &rValues, c10::optional< at::Tensor > out)
 Resample the given 3D sinogram at locations corresponding to the given 2D sinogram grid (phiValues, rValues), according to the 2D-3D image registration method based on Grangeat's relation.
 
__host__ at::Tensor reg23::ResampleSinogram3D_CUDA (const at::Tensor &sinogram3d, const std::string &sinogramType, double rSpacing, const at::Tensor &projectionMatrix, const at::Tensor &phiValues, const at::Tensor &rValues, c10::optional< at::Tensor > out)
 An implementation of reg23::ResampleSinogram3D_CPU that uses CUDA parallelisation.
 
__host__ at::Tensor reg23::ResampleSinogram3DCUDATexture (int64_t sinogram3dTextureHandle, int64_t sinogramWidth, int64_t sinogramHeight, int64_t sinogramDepth, const std::string &sinogramType, double rSpacing, const at::Tensor &projectionMatrix, const at::Tensor &phiValues, const at::Tensor &rValues, c10::optional< at::Tensor > out)
 An implementation of reg23::ResampleSinogram3D_CUDA that takes a handle to a pre-allocated CUDA texture instead of a PyTorch tensor.
 
std::tuple< at::Tensor, double, double, double, double, doublereg23::NormalisedCrossCorrelation_CPU (const at::Tensor &a, const at::Tensor &b)
 Additionally returns intermediate quantities useful for evaluating the backward pass.
 
__host__ std::tuple< at::Tensor, double, double, double, double, doublereg23::NormalisedCrossCorrelation_CUDA (const at::Tensor &a, const at::Tensor &b)
 An implementation of reg23::NormalisedCrossCorrelation_CPU that uses CUDA parallelisation.
 

Detailed Description

Functions that have Python bindings generated and are accessible in the PyTorch extension.

Function Documentation

◆ DRadon2DDR_CPU()

at::Tensor reg23::DRadon2DDR_CPU ( const at::Tensor &  image,
const at::Tensor &  imageSpacing,
const at::Tensor &  phiValues,
const at::Tensor &  rValues,
int64_t  samplesPerLine 
)

Compute the derivative with respect to plane-origin distance of an approximation of the Radon transform of the given 2D image.

Parameters
imagea tensor of size (N, M) containing torch.float32s: The input image to take the Radon transform of
imageSpacinga tensor of size (2,): The spacing between the image columns and rows
phiValuesa tensor of any size containing torch.float32s: The values of the polar coordinate phi at which to calculate approximations for line integrals through the given image
rValuesa tensor of the same size as phiValues containing torch.float32s: The values of the polar coordinate r at which to calculate approximations for line integrals through the given image
samplesPerLineThe number of samples to take from the image for approximating each line integral
Returns
a tensor matching phiValues in size containing torch.float32s: The derivatives with respect to plane-origin distance of the approximations (according to the Radon2D_... functions) of the line integrals through the given image at the given polar coordinate locations

The polar coordinates should be defined according to the convention stated in reg23/Conventions.md

◆ DRadon2DDR_CUDA()

__host__ at::Tensor reg23::DRadon2DDR_CUDA ( const at::Tensor &  image,
const at::Tensor &  imageSpacing,
const at::Tensor &  phiValues,
const at::Tensor &  rValues,
int64_t  samplesPerLine 
)

An implementation of reg23::DRadon2DDR_CPU that uses CUDA parallelisation.

A single kernel launch is made, with each kernel calculating one line integral approximation.

◆ DRadon3DDR_CPU()

at::Tensor reg23::DRadon3DDR_CPU ( const at::Tensor &  volume,
const at::Tensor &  volumeSpacing,
const at::Tensor &  phiValues,
const at::Tensor &  thetaValues,
const at::Tensor &  rValues,
int64_t  samplesPerDirection 
)

Compute the derivative with respect to plane-origin distance of an approximation of the Radon transform of the given 3D volume.

Parameters
volumea tensor of size (P,Q,R) containing torch.float32s: The input volume to take the Radon transform of
volumeSpacinga tensor of size (3,): The spacing between the volume layers in each cartesian direction
phiValuesa tensor of any size containing torch.float32s: The values of the spherical coordinate phi at which to calculate the derivatives of approximations for plane integrals through the given volume
thetaValuesa tensor of the same size as phiValues containing torch.float32s: The values of the spherical coordinate theta at which to calculate the derivatives of approximations for plane integrals through the given volume
rValuesa tensor of the same size as phiValues containing torch.float32s: The values of the spherical coordinate r at which to calculate the derivatives of approximations for plane integrals through the given volume
samplesPerDirectionThe square-root of the number of samples to take from the volume for approximating each plane integral
Returns
a tensor of the same size as phiValues containing torch.float32s: The derivatives with respect to plane-origin distance of the approximations (according to the Radon3D_... functions) of the plane integrals through the given volume at the given spherical coordinate locations

The spherical coordinates should be defined according to the convention stated in reg23/Conventions.md

◆ DRadon3DDR_CUDA()

__host__ at::Tensor reg23::DRadon3DDR_CUDA ( const at::Tensor &  volume,
const at::Tensor &  volumeSpacing,
const at::Tensor &  phiValues,
const at::Tensor &  thetaValues,
const at::Tensor &  rValues,
int64_t  samplesPerDirection 
)

An implementation of reg23::DRadon3DDR_CPU that uses CUDA parallelisation.

A single kernel launch is made, with each kernel calculating one plane integral approximation.

◆ GridSample3D_CPU()

at::Tensor reg23::GridSample3D_CPU ( const at::Tensor &  input,
const at::Tensor &  grid,
const std::string &  addressModeX,
const std::string &  addressModeY,
const std::string &  addressModeZ,
c10::optional< at::Tensor >  out 
)

Sample the given 3D input tensor at the positions given in grid according to the given address mode using bilinear interpolation. This implementation is single-threaded.

Parameters
inputa tensor of size (P,Q,R) containing torch.float32s: The volume to sample
grida tensor of size (..., 3) containing torch.float32s: The grid to sample at; the last dimension represents 3D locations within the input volume between (-1, -1, -1) and (1, 1, 1)
addressModeXeither "wrap" or "zero": The address mode used for dimension X
addressModeYeither "wrap" or "zero": The address mode used for dimension Y
addressModeZeither "wrap" or "zero": The address mode used for dimension Z
Returns
a tensor of size (...) containing torch.float32s: The sampled values in the given volume at the locations in the given grid; this will be the same size as grid, minus the last dimension.

Address modes

  • "zero": sampling locations outside (-1, -1, -1) and (1, 1, 1) will be read as 0
  • "wrap": sampling locations (x, y, z) outside (-1, -1, -1) and (1, 1, 1) will be read wrapped back according to ((x + 1) mod 2 - 1, (y + 1) mod 2 - 1, (z + 1) mod 2 - 1)

◆ GridSample3D_CUDA()

__host__ at::Tensor reg23::GridSample3D_CUDA ( const at::Tensor &  input,
const at::Tensor &  grid,
const std::string &  addressModeX,
const std::string &  addressModeY,
const std::string &  addressModeZ,
c10::optional< at::Tensor >  out 
)

An implementation of reg23::GridSample3D_CPU that uses CUDA parallelisation.

◆ NormalisedCrossCorrelation_CPU()

std::tuple< at::Tensor, double, double, double, double, double > reg23::NormalisedCrossCorrelation_CPU ( const at::Tensor &  a,
const at::Tensor &  b 
)

Additionally returns intermediate quantities useful for evaluating the backward pass.

Parameters
aa tensor of any size containing torch.float32s
ba tensor of the same size as a containing torch.float32s
Returns
a tuple: (ZNCC; sum of values in a; sum of values in b; the numerator of the fraction shown below; the left sqrt in the denominator of the fraction shown below; the right sqrt in the denominator of the fraction shown below)

ZNCC = (N * sum(a_i * b_i) - sum(a_i) * sum(b_i)) / (sqrt(N * sum(a_i^2) - sum(a_i)^2) * sqrt(N * sum(b_i^2) - sum(b_i)^2))

◆ NormalisedCrossCorrelation_CUDA()

__host__ std::tuple< at::Tensor, double, double, double, double, double > reg23::NormalisedCrossCorrelation_CUDA ( const at::Tensor &  a,
const at::Tensor &  b 
)

An implementation of reg23::NormalisedCrossCorrelation_CPU that uses CUDA parallelisation.

◆ ProjectDRR_backward_CPU()

at::Tensor reg23::ProjectDRR_backward_CPU ( const at::Tensor &  volume,
const at::Tensor &  voxelSpacing,
const at::Tensor &  homographyMatrixInverse,
double  sourceDistance,
int64_t  outputWidth,
int64_t  outputHeight,
const at::Tensor &  outputOffset,
const at::Tensor &  detectorSpacing,
const at::Tensor &  dLossDDRR 
)

Evaluate the derivative of some scalar loss that is a function of a DRR projected from the given volume at the given transformation, with respect to the inverse homography matrix.

Parameters
volumea tensor of size (P,Q,R): The CT volume through which to project the DRR
voxelSpacinga tensor of size (3,): The spacing in [mm] between the volume layers in each cartesian direction
homographyMatrixInversea tensor of size (4, 4): The column-major matrix representing the homography transformation of the volume.
sourceDistanceThe distance in [mm] of the source from the detector array
outputWidthThe width in pixels of the DRR to generate.
outputHeightThe height in pixels of the DRR to generate.
outputOffseta tensor of size (2,) containing torch.float64s: The offset in [mm] of the centre of the DRR from the central ray from the X-ray source perpendicularly onto the detector array.
detectorSpacinga tensor of size (2,): The spacing in [mm] between the columns and rows of the DRR image.
dLossDDRRa tensor of size (outputHeight, outputWidth): The derivative of the loss w.r.t. the projected DRR image
Returns
tensor of size (4, 4): The derivative of the loss w.r.t. the inverse homography matrix

◆ ProjectDRR_backward_CUDA()

__host__ at::Tensor reg23::ProjectDRR_backward_CUDA ( const at::Tensor &  volume,
const at::Tensor &  voxelSpacing,
const at::Tensor &  homographyMatrixInverse,
double  sourceDistance,
int64_t  outputWidth,
int64_t  outputHeight,
const at::Tensor &  outputOffset,
const at::Tensor &  detectorSpacing,
const at::Tensor &  dLossDDRR 
)

An implementation of reg23::ProjectDRR_backward_CPU that uses CUDA parallelisation.

◆ ProjectDRR_CPU()

at::Tensor reg23::ProjectDRR_CPU ( const at::Tensor &  volume,
const at::Tensor &  voxelSpacing,
const at::Tensor &  homographyMatrixInverse,
double  sourceDistance,
int64_t  outputWidth,
int64_t  outputHeight,
const at::Tensor &  outputOffset,
const at::Tensor &  detectorSpacing 
)

Generate a DRR from the given volume at the given transformation.

Parameters
volumea tensor of size (P,Q,R): The CT volume through which to project the DRR
voxelSpacinga tensor of size (3,): The spacing in [mm] between the volume layers in each cartesian direction
homographyMatrixInversea tensor of size (4, 4): The column-major matrix representing the homography transformation of the volume.
sourceDistanceThe distance in [mm] of the source from the detector array
outputWidthThe width in pixels of the DRR to generate.
outputHeightThe height in pixels of the DRR to generate.
outputOffseta tensor of size (2,) containing torch.float64s: The offset in [mm] of the centre of the DRR from the central ray from the X-ray source perpendicularly onto the detector array.
detectorSpacinga tensor of size (2,): The spacing in [mm] between the columns and rows of the DRR image.
Returns
a tensor of size (outputHeight, outputWidth): The DRR projection through the given volume at the given transformation.

The number of samples taken along each ray through the volume for approximation of each ray integral is equal to the largest dimension of the given volume.

◆ ProjectDRR_CUDA()

__host__ at::Tensor reg23::ProjectDRR_CUDA ( const at::Tensor &  volume,
const at::Tensor &  voxelSpacing,
const at::Tensor &  homographyMatrixInverse,
double  sourceDistance,
int64_t  outputWidth,
int64_t  outputHeight,
const at::Tensor &  outputOffset,
const at::Tensor &  detectorSpacing 
)

An implementation of reg23::ProjectDRR_CPU that uses CUDA parallelisation.

◆ ProjectDRRCuboidMask_CPU()

at::Tensor reg23::ProjectDRRCuboidMask_CPU ( const at::Tensor &  volumeSize,
const at::Tensor &  voxelSpacing,
const at::Tensor &  homographyMatrixInverse,
double  sourceDistance,
int64_t  outputWidth,
int64_t  outputHeight,
const at::Tensor &  outputOffset,
const at::Tensor &  detectorSpacing 
)

Calculate the distances of the intersections between the DRR-generation rays and the domain of the volume data.

Parameters
volumeSizea tensor of size (3,): the sizes of the volume data tensor in the order (X, Y, Z) - this is the reverse order of tensor.size()
voxelSpacinga tensor of size (3,): the spacing in the (X, Y, Z) directions between adjactent voxels in the volume
homographyMatrixInversea tensor of size (4, 4): the inverse homography matrix of the volume transformation
sourceDistancethe distance of the X-ray source from the detector array
outputWidththe width of the detector array
outputHeightthe height of the detector array
outputOffseta tensor of size (2,) containing torch.float64s: The offset in [mm] of the centre of the DRR from the central ray from the X-ray source perpendicularly onto the detector array.
detectorSpacinga tensor of size (2,): The spacing in [mm] between the columns and rows of the DRR image.
Returns
A tensor of size (outputHeight, outputWidth): The distance each DRR ray travels through the volume.

◆ ProjectDRRCuboidMask_CUDA()

__host__ at::Tensor reg23::ProjectDRRCuboidMask_CUDA ( const at::Tensor &  volumeSize,
const at::Tensor &  voxelSpacing,
const at::Tensor &  homographyMatrixInverse,
double  sourceDistance,
int64_t  outputWidth,
int64_t  outputHeight,
const at::Tensor &  outputOffset,
const at::Tensor &  detectorSpacing 
)

An implementation of reg23::ProjectDRRCuboidMask_CPU that uses CUDA parallelisation.

◆ ProjectDRRsBatched_CUDA()

__host__ at::Tensor reg23::ProjectDRRsBatched_CUDA ( const at::Tensor &  volume,
const at::Tensor &  voxelSpacing,
const at::Tensor &  invHMatrices,
double  sourceDistance,
int64_t  outputWidth,
int64_t  outputHeight,
const at::Tensor &  outputOffset,
const at::Tensor &  detectorSpacing 
)

An implementation similar to reg23::ProjectDRR_CUDA that evaluates projections for multiple transformations in parallel.

Parameters
volumea tensor of size (P,Q,R): The CT volume through which to project the DRR
voxelSpacinga tensor of size (3,): The spacing in [mm] between the volume layers in each cartesian direction
invHMatricesa tensor of size (..., 4, 4): The column-major matrices representing the homography transformations of the volume.
sourceDistanceThe distance in [mm] of the source from the detector array
outputWidthThe width in pixels of the DRR to generate.
outputHeightThe height in pixels of the DRR to generate.
outputOffseta tensor of size (2,) containing torch.float64s: The offset in [mm] of the centre of the DRR from the central ray from the X-ray source perpendicularly onto the detector array.
detectorSpacinga tensor of size (2,): The spacing in [mm] between the columns and rows of the DRR image.
Returns
a tensor of size (..., outputHeight, outputWidth): The DRR projection through the given volume at the given transformation. The '...' dimensions correspond to those of the invHMatrices tensor.

◆ Radon2D_CPU()

at::Tensor reg23::Radon2D_CPU ( const at::Tensor &  image,
const at::Tensor &  imageSpacing,
const at::Tensor &  phiValues,
const at::Tensor &  rValues,
int64_t  samplesPerLine 
)

Compute an approximation of the Radon transform of the given 2D image.

Parameters
imagea tensor of size (N, M) containing torch.float32s: The input image to take the Radon transform of
imageSpacinga tensor of size (2,): The spacing between the image columns and rows
phiValuesa tensor of any size containing torch.float32s: The values of the polar coordinate phi at which to calculate approximations for line integrals through the given image
rValuesa tensor of the same size as phiValues containing torch.float32s: The values of the polar coordinate r at which to calculate approximations for line integrals through the given image
samplesPerLineThe number of samples to take from the image for approximating each line integral
Returns
a tensor matching phiValues in size containing torch.float32s: Approximations of the line integrals through the given image at the given polar coordinate locations

The polar coordinates should be defined according to the convention stated in reg23/Conventions.md

◆ Radon2D_CUDA()

__host__ at::Tensor reg23::Radon2D_CUDA ( const at::Tensor &  image,
const at::Tensor &  imageSpacing,
const at::Tensor &  phiValues,
const at::Tensor &  rValues,
int64_t  samplesPerLine 
)

An implementation of reg23::Radon2D_CPU that uses CUDA parallelisation.

A single kernel launch is made, with each kernel calculating one line integral approximation.

◆ Radon3D_CPU()

at::Tensor reg23::Radon3D_CPU ( const at::Tensor &  volume,
const at::Tensor &  volumeSpacing,
const at::Tensor &  phiValues,
const at::Tensor &  thetaValues,
const at::Tensor &  rValues,
int64_t  samplesPerDirection 
)

Compute an approximation of the Radon transform of the given 3D volume.

Parameters
volumea tensor of size (P,Q,R) containing torch.float32s: The input volume to take the Radon transform of
volumeSpacinga tensor of size (3,): The spacing between the volume layers in each cartesian direction
phiValuesa tensor of any size containing torch.float32s: The values of the spherical coordinate $phi$ at which to calculate approximations for plane integrals through the given volume
thetaValuesa tensor of the same size as phiValues containing torch.float32s: The values of the spherical coordinate $theta$ at which to calculate approximations for plane integrals through the given volume
rValuesa tensor of the same size as phiValues containing torch.float32s: The values of the spherical coordinate $r$ at which to calculate approximations for plane integrals through the given volume
samplesPerDirectionThe square-root of the number of samples to take from the volume for approximating each plane integral
Returns
a tensor of the same size as phiValues containing torch.float32s: Approximations of the plane integrals through the given volume at the given spherical coordinate locations

The spherical coordinates should be defined according to the convention stated in reg23/Conventions.md

◆ Radon3D_CUDA()

__host__ at::Tensor reg23::Radon3D_CUDA ( const at::Tensor &  volume,
const at::Tensor &  volumeSpacing,
const at::Tensor &  phiValues,
const at::Tensor &  thetaValues,
const at::Tensor &  rValues,
int64_t  samplesPerDirection 
)

An implementation of reg23::Radon3D_CPU that uses CUDA parallelisation.

A single kernel launch is made, with each kernel calculating one plane integral approximation.

◆ ResampleSinogram3D_CPU()

at::Tensor reg23::ResampleSinogram3D_CPU ( const at::Tensor &  sinogram3d,
const std::string &  sinogramType,
double  rSpacing,
const at::Tensor &  projectionMatrix,
const at::Tensor &  phiValues,
const at::Tensor &  rValues,
c10::optional< at::Tensor >  out 
)

Resample the given 3D sinogram at locations corresponding to the given 2D sinogram grid (phiValues, rValues), according to the 2D-3D image registration method based on Grangeat's relation.

Parameters
sinogram3da tensor of size (P,Q,R) containing torch.float32s: The 3D sinogram to resample. It must contain plane integrals evenly spaced along phi, theta and r over each dimension respectively. This is the pre-calculated derivative of the Radon transform of the CT volume.
sinogramTypea string; one of the following: "classic", "healpix". see documentation for enum reg23::ResampleSinogram3D::SinogramType
rSpacingThe spacing in r between the plane integrals in the sinogram.
projectionMatrixa tensor of size (4, 4) containing torch.float32s: The projection matrix PH at the desired transformation.
phiValuesa tensor of any size containing torch.float32s: The polar coordinate phi values of the fixed image at which to resample corresponding values in the given sinogram
rValuesa tensor of the same size as phiValues containing torch.float32s: The polar coordinate r values of the fixed image at which to resample corresponding values in the given sinogram
Returns
a tensor of the same size as phiValues containing torch.float32s: The values sampled from the given volume sinogram at the locations corresponding to the given polar coordinates (phiValues, rValues), according to the 2D-3D image registration method based on Grangeat's relation:

Frysch R, Pfeiffer T, Rose G. A novel approach to 2D/3D registration of X-ray images using Grangeat's relation. Med Image Anal. 2021 Jan;67:101815. doi: 10.1016/j.media.2020.101815. Epub 2020 Sep 30. PMID: 33065470.

Note: Assumes that the projection matrix projects onto the x-y plane, and that the radial coordinates (phi, r) in that plane measure phi right-hand rule about the z-axis from the positive x-direction

◆ ResampleSinogram3D_CUDA()

__host__ at::Tensor reg23::ResampleSinogram3D_CUDA ( const at::Tensor &  sinogram3d,
const std::string &  sinogramType,
double  rSpacing,
const at::Tensor &  projectionMatrix,
const at::Tensor &  phiValues,
const at::Tensor &  rValues,
c10::optional< at::Tensor >  out 
)

An implementation of reg23::ResampleSinogram3D_CPU that uses CUDA parallelisation.

◆ ResampleSinogram3DCUDATexture()

__host__ at::Tensor reg23::ResampleSinogram3DCUDATexture ( int64_t  sinogram3dTextureHandle,
int64_t  sinogramWidth,
int64_t  sinogramHeight,
int64_t  sinogramDepth,
const std::string &  sinogramType,
double  rSpacing,
const at::Tensor &  projectionMatrix,
const at::Tensor &  phiValues,
const at::Tensor &  rValues,
c10::optional< at::Tensor >  out 
)

An implementation of reg23::ResampleSinogram3D_CUDA that takes a handle to a pre-allocated CUDA texture instead of a PyTorch tensor.