Skip to content
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
26 changes: 26 additions & 0 deletions applications/rtkconjugategradient/rtkconjugategradient.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -29,8 +29,10 @@

#ifdef RTK_USE_CUDA
# include <itkCudaImage.h>
# include <itkCovariantVector.h>
#endif
#include <itkImageFileWriter.h>
#include <itkImageFileReader.h>

int
main(int argc, char * argv[])
Expand Down Expand Up @@ -109,6 +111,30 @@ main(int argc, char * argv[])
auto conjugategradient = rtk::ConjugateGradientConeBeamReconstructionFilter<OutputImageType>::New();
SetForwardProjectionFromGgo(args_info, conjugategradient.GetPointer());
SetBackProjectionFromGgo(args_info, conjugategradient.GetPointer());


#ifdef RTK_USE_CUDA
{
const bool warpRequested = (args_info.fp_arg == fp_arg_CudaWarpRayCast) || (args_info.bp_arg == bp_arg_CudaWarp);
if (warpRequested && !args_info.dvf_given)
{
itkGenericExceptionMacro(<< "Warp projectors selected but no --dvf provided. Please provide "
"a DVF in physical units (mm).");
}
if (args_info.dvf_given)
{
using VectorType = itk::CovariantVector<OutputPixelType, 3>;
using DVF3DCudaType = itk::CudaImage<VectorType, 3>;

auto dvf3dReader = itk::ImageFileReader<DVF3DCudaType>::New();
dvf3dReader->SetFileName(args_info.dvf_arg);
TRY_AND_EXIT_ON_ITK_EXCEPTION(dvf3dReader->Update())

conjugategradient->SetDisplacementField(reinterpret_cast<const OutputImageType *>(dvf3dReader->GetOutput()));
}
}
#endif

conjugategradient->SetInputVolume(inputFilter->GetOutput());
conjugategradient->SetInputProjectionStack(reader->GetOutput());
conjugategradient->SetInputWeights(inputWeights);
Expand Down
1 change: 1 addition & 0 deletions applications/rtkconjugategradient/rtkconjugategradient.ggo
Original file line number Diff line number Diff line change
Expand Up @@ -13,3 +13,4 @@ option "tikhonov" - "Tikhonov regularization weight"
option "nocudacg" - "Do not perform conjugate gradient calculations on GPU" flag off
option "mask" m "Apply a support binary mask: reconstruction kept null outside the mask)" string no
option "nodisplaced" - "Disable the displaced detector filter" flag off
option "dvf" - "3D DVF" string no
2 changes: 1 addition & 1 deletion applications/rtkiterativefdk/rtkiterativefdk.ggo
Original file line number Diff line number Diff line change
Expand Up @@ -16,7 +16,7 @@ option "hann" - "Cut frequency for hann window in ]0,1] (0.0 disables it)"
option "hannY" - "Cut frequency for hann window in ]0,1] (0.0 disables it)" double no default="0.0"

section "Projectors"
option "fp" f "Forward projection method" values="Joseph","CudaRayCast","JosephAttenuated","Zeng" enum no default="Joseph"
option "fp" f "Forward projection method" values="Joseph","CudaRayCast","JosephAttenuated","Zeng","CudaWarpRayCast" enum no default="Joseph"
option "attenuationmap" - "Attenuation map relative to the volume to perfom the attenuation correction" string no
option "sigmazero" - "PSF value at a distance of 0 meter of the detector" double no
option "alphapsf" - "Slope of the PSF against the detector distance" double no
Expand Down
6 changes: 3 additions & 3 deletions applications/rtkprojectors_section.ggo
Original file line number Diff line number Diff line change
@@ -1,7 +1,7 @@
section "Projectors"
option "fp" f "Forward projection method" values="Joseph","CudaRayCast","JosephAttenuated","Zeng" enum no default="Joseph"
option "bp" b "Back projection method" values="VoxelBasedBackProjection","Joseph","CudaVoxelBased","CudaRayCast","JosephAttenuated", "Zeng" enum no default="VoxelBasedBackProjection"
option "step" - "Step size along ray (for CudaRayCast only)" double no default="1"
option "fp" f "Forward projection method" values="Joseph","CudaRayCast","JosephAttenuated","Zeng","CudaWarpRayCast" enum no default="Joseph"
option "bp" b "Back projection method" values="VoxelBasedBackProjection","Joseph","CudaVoxelBased","CudaRayCast","JosephAttenuated","Zeng","CudaWarp" enum no default="VoxelBasedBackProjection"
option "step" - "Step size along ray (for CudaRayCast and CudaWarpRayCast)" double no default="1"
option "attenuationmap" - "Attenuation map relative to the volume to perfom the attenuation correction (JosephAttenuated and Zeng)" string no
option "sigmazero" - "PSF value at a distance of 0 meter of the detector (Zeng only)" double no
option "alphapsf" - "Slope of the PSF against the detector distance (Zeng only)" double no
Expand Down
8 changes: 8 additions & 0 deletions include/rtkGgoFunctions.h
Original file line number Diff line number Diff line change
Expand Up @@ -331,6 +331,9 @@ SetBackProjectionFromGgo(const TArgsInfo & args_info, TIterativeReconstructionFi
if (args_info.attenuationmap_given)
recon->SetAttenuationMap(attenuationMap);
break;
case (6): // bp_arg_CudaWarpBackProjection
recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_CUDAWARP);
break;
}
}

Expand Down Expand Up @@ -400,6 +403,11 @@ SetForwardProjectionFromGgo(const TArgsInfo & args_info, TIterativeReconstructio
if (args_info.attenuationmap_given)
recon->SetAttenuationMap(attenuationMap);
break;
case (4): // fp_arg_CudaWarpRayCast
recon->SetForwardProjectionFilter(TIterativeReconstructionFilter::FP_CUDAWARPRAYCAST);
if (args_info.step_given)
recon->SetStepSize(args_info.step_arg);
break;
}
}

Expand Down
85 changes: 83 additions & 2 deletions include/rtkIterativeConeBeamReconstructionFilter.h
Original file line number Diff line number Diff line change
Expand Up @@ -33,6 +33,8 @@
# include "rtkCudaForwardProjectionImageFilter.h"
# include "rtkCudaBackProjectionImageFilter.h"
# include "rtkCudaRayCastBackProjectionImageFilter.h"
# include "rtkCudaWarpForwardProjectionImageFilter.h"
# include "rtkCudaWarpBackProjectionImageFilter.h"
#endif

#include <random>
Expand Down Expand Up @@ -74,16 +76,19 @@ class ITK_TEMPLATE_EXPORT IterativeConeBeamReconstructionFilter
FP_JOSEPH = 0,
FP_CUDARAYCAST = 2,
FP_JOSEPHATTENUATED = 3,
FP_ZENG = 4
FP_ZENG = 4,
FP_CUDAWARPRAYCAST = 5
} ForwardProjectionType;
typedef enum
{
BP_VOXELBASED = 0,
BP_JOSEPH = 1,
BP_CUDAVOXELBASED = 2,
BP_CUDAWARP = 3,
BP_CUDARAYCAST = 4,
BP_JOSEPHATTENUATED = 5,
BP_ZENG = 6
BP_ZENG = 6,

} BackProjectionType;

/** Typedefs of each subfilter of this composite filter */
Expand Down Expand Up @@ -158,6 +163,18 @@ class ITK_TEMPLATE_EXPORT IterativeConeBeamReconstructionFilter
return static_cast<const TClipImageType *>(this->itk::ProcessObject::GetInput("SuperiorClipImage"));
}

/** Set/Get the displacement field (DVF) required by warp forward/back projectors. */
void
SetDisplacementField(const VolumeType * dvf)
{
this->SetInput("DisplacementField", const_cast<VolumeType *>(dvf));
}
typename VolumeType::ConstPointer
GetDisplacementField()
{
return static_cast<const VolumeType *>(this->itk::ProcessObject::GetInput("DisplacementField"));
}

/** Get / Set the sigma zero of the PSF. Default is 1.5417233052142099 */
itkGetMacro(SigmaZero, double);
itkSetMacro(SigmaZero, double);
Expand Down Expand Up @@ -367,6 +384,70 @@ class ITK_TEMPLATE_EXPORT IterativeConeBeamReconstructionFilter
}


template <typename ImageType, EnableCudaScalarType<ImageType> * = nullptr>
ForwardProjectionPointerType
InstantiateCudaWarpForwardProjection()
{
ForwardProjectionPointerType fw;
#ifdef RTK_USE_CUDA
fw = CudaWarpForwardProjectionImageFilter::New();
dynamic_cast<rtk::CudaWarpForwardProjectionImageFilter *>(fw.GetPointer())->SetStepSize(m_StepSize);
if (this->GetDisplacementField().IsNull())
{
itkExceptionMacro(
<< "CudaWarpForwardProjectionImageFilter requires a displacement field. Call SetDisplacementField().");
}
else
{
dynamic_cast<rtk::CudaWarpForwardProjectionImageFilter *>(fw.GetPointer())
->SetDisplacementField(reinterpret_cast<typename rtk::CudaWarpForwardProjectionImageFilter::DVFType *>(
const_cast<VolumeType *>(this->GetDisplacementField().GetPointer())));
}
#endif
return fw;
}


template <typename ImageType, DisableCudaScalarType<ImageType> * = nullptr>
ForwardProjectionPointerType
InstantiateCudaWarpForwardProjection()
{
itkGenericExceptionMacro(<< "CudaWarpForwardProjectionImageFilter only available with 3D CudaImage of float.");
return nullptr;
}


template <typename ImageType, EnableCudaScalarType<ImageType> * = nullptr>
BackProjectionPointerType
InstantiateCudaWarpBackProjection()
{
BackProjectionPointerType bp;
#ifdef RTK_USE_CUDA
bp = CudaWarpBackProjectionImageFilter::New();
if (this->GetDisplacementField().IsNull())
{
itkExceptionMacro(
<< "CudaWarpBackProjectionImageFilter requires a displacement field. Call SetDisplacementField().");
}
else
{
dynamic_cast<rtk::CudaWarpBackProjectionImageFilter *>(bp.GetPointer())
->SetDisplacementField(reinterpret_cast<typename rtk::CudaWarpBackProjectionImageFilter::DVFType *>(
const_cast<VolumeType *>(this->GetDisplacementField().GetPointer())));
}
#endif
return bp;
}

template <typename ImageType, DisableCudaScalarType<ImageType> * = nullptr>
BackProjectionPointerType
InstantiateCudaWarpBackProjection()
{
itkGenericExceptionMacro(<< "CudaWarpBackProjectionImageFilter only available with 3D CudaImage of float.");
return nullptr;
}


template <typename ImageType, EnableVectorType<ImageType> * = nullptr>
BackProjectionPointerType
InstantiateJosephBackAttenuatedProjection()
Expand Down
6 changes: 6 additions & 0 deletions include/rtkIterativeConeBeamReconstructionFilter.hxx
Original file line number Diff line number Diff line change
Expand Up @@ -59,6 +59,9 @@ IterativeConeBeamReconstructionFilter<TOutputImage, ProjectionStackType>::Instan
case (FP_ZENG):
fw = InstantiateZengForwardProjection<ProjectionStackType>();
break;
case (FP_CUDAWARPRAYCAST):
fw = InstantiateCudaWarpForwardProjection<ProjectionStackType>();
break;
default:
itkGenericExceptionMacro(<< "Unhandled --fp value.");
}
Expand Down Expand Up @@ -90,6 +93,9 @@ IterativeConeBeamReconstructionFilter<TOutputImage, ProjectionStackType>::Instan
case (BP_ZENG):
bp = InstantiateZengBackProjection<ProjectionStackType>();
break;
case (BP_CUDAWARP):
bp = InstantiateCudaWarpBackProjection<ProjectionStackType>();
break;
default:
itkGenericExceptionMacro(<< "Unhandled --bp value.");
}
Expand Down
21 changes: 21 additions & 0 deletions test/rtkconjugategradientreconstructiontest.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -7,6 +7,7 @@
#ifdef USE_CUDA
# include "itkCudaImage.h"
#endif
#include <itkCovariantVector.h>

/**
* \file rtkconjugategradientreconstructiontest.cxx
Expand Down Expand Up @@ -139,6 +140,26 @@ main(int, char **)

CheckImageQuality<OutputImageType>(conjugategradient->GetOutput(), dsl->GetOutput(), 0.08, 23, 2.0);
std::cout << "\n\nTest PASSED! " << std::endl;

std::cout << "\n\n****** Case 3b: CUDA Warp Forward/Back Projectors ******" << std::endl;
using DVFPixelType = itk::CovariantVector<OutputPixelType, 3>;
using DVFImageType = itk::CudaImage<DVFPixelType, 3>;
auto dvfSource = rtk::ConstantImageSource<DVFImageType>::New();
dvfSource->SetOrigin(tomographySource->GetOutput()->GetOrigin());
dvfSource->SetSpacing(tomographySource->GetOutput()->GetSpacing());
dvfSource->SetSize(tomographySource->GetOutput()->GetLargestPossibleRegion().GetSize());
DVFPixelType dvf;
dvf.Fill(0.0);
dvfSource->SetConstant(dvf);
TRY_AND_EXIT_ON_ITK_EXCEPTION(dvfSource->Update())

conjugategradient->SetDisplacementField(reinterpret_cast<const OutputImageType *>(dvfSource->GetOutput()));
conjugategradient->SetForwardProjectionFilter(ConjugateGradientType::FP_CUDAWARPRAYCAST);
conjugategradient->SetBackProjectionFilter(ConjugateGradientType::BP_CUDAWARP);
TRY_AND_EXIT_ON_ITK_EXCEPTION(conjugategradient->Update());

CheckImageQuality<OutputImageType>(conjugategradient->GetOutput(), dsl->GetOutput(), 0.08, 23, 2.0);
std::cout << "\n\nTest PASSED! " << std::endl;
#endif

std::cout << "\n\n****** Case 4: Joseph Backprojector, weighted least squares ******" << std::endl;
Expand Down
Loading