diff --git a/applications/rtkconjugategradient/rtkconjugategradient.cxx b/applications/rtkconjugategradient/rtkconjugategradient.cxx index 6eec4ee37..b1e7c64ca 100644 --- a/applications/rtkconjugategradient/rtkconjugategradient.cxx +++ b/applications/rtkconjugategradient/rtkconjugategradient.cxx @@ -29,8 +29,10 @@ #ifdef RTK_USE_CUDA # include +# include #endif #include +#include int main(int argc, char * argv[]) @@ -109,6 +111,30 @@ main(int argc, char * argv[]) auto conjugategradient = rtk::ConjugateGradientConeBeamReconstructionFilter::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; + using DVF3DCudaType = itk::CudaImage; + + auto dvf3dReader = itk::ImageFileReader::New(); + dvf3dReader->SetFileName(args_info.dvf_arg); + TRY_AND_EXIT_ON_ITK_EXCEPTION(dvf3dReader->Update()) + + conjugategradient->SetDisplacementField(reinterpret_cast(dvf3dReader->GetOutput())); + } + } +#endif + conjugategradient->SetInputVolume(inputFilter->GetOutput()); conjugategradient->SetInputProjectionStack(reader->GetOutput()); conjugategradient->SetInputWeights(inputWeights); diff --git a/applications/rtkconjugategradient/rtkconjugategradient.ggo b/applications/rtkconjugategradient/rtkconjugategradient.ggo index 864b490a4..bba2ec897 100644 --- a/applications/rtkconjugategradient/rtkconjugategradient.ggo +++ b/applications/rtkconjugategradient/rtkconjugategradient.ggo @@ -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 \ No newline at end of file diff --git a/applications/rtkiterativefdk/rtkiterativefdk.ggo b/applications/rtkiterativefdk/rtkiterativefdk.ggo index f6e8ff868..1b8ae0e7a 100644 --- a/applications/rtkiterativefdk/rtkiterativefdk.ggo +++ b/applications/rtkiterativefdk/rtkiterativefdk.ggo @@ -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 diff --git a/applications/rtkprojectors_section.ggo b/applications/rtkprojectors_section.ggo index 1f33fb0ea..acbdf6cc3 100644 --- a/applications/rtkprojectors_section.ggo +++ b/applications/rtkprojectors_section.ggo @@ -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 diff --git a/include/rtkGgoFunctions.h b/include/rtkGgoFunctions.h index 79d07117a..34feaa556 100644 --- a/include/rtkGgoFunctions.h +++ b/include/rtkGgoFunctions.h @@ -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; } } @@ -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; } } diff --git a/include/rtkIterativeConeBeamReconstructionFilter.h b/include/rtkIterativeConeBeamReconstructionFilter.h index 589566acb..1a05a17b9 100644 --- a/include/rtkIterativeConeBeamReconstructionFilter.h +++ b/include/rtkIterativeConeBeamReconstructionFilter.h @@ -33,6 +33,8 @@ # include "rtkCudaForwardProjectionImageFilter.h" # include "rtkCudaBackProjectionImageFilter.h" # include "rtkCudaRayCastBackProjectionImageFilter.h" +# include "rtkCudaWarpForwardProjectionImageFilter.h" +# include "rtkCudaWarpBackProjectionImageFilter.h" #endif #include @@ -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 */ @@ -158,6 +163,18 @@ class ITK_TEMPLATE_EXPORT IterativeConeBeamReconstructionFilter return static_cast(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(dvf)); + } + typename VolumeType::ConstPointer + GetDisplacementField() + { + return static_cast(this->itk::ProcessObject::GetInput("DisplacementField")); + } + /** Get / Set the sigma zero of the PSF. Default is 1.5417233052142099 */ itkGetMacro(SigmaZero, double); itkSetMacro(SigmaZero, double); @@ -367,6 +384,70 @@ class ITK_TEMPLATE_EXPORT IterativeConeBeamReconstructionFilter } + template * = nullptr> + ForwardProjectionPointerType + InstantiateCudaWarpForwardProjection() + { + ForwardProjectionPointerType fw; +#ifdef RTK_USE_CUDA + fw = CudaWarpForwardProjectionImageFilter::New(); + dynamic_cast(fw.GetPointer())->SetStepSize(m_StepSize); + if (this->GetDisplacementField().IsNull()) + { + itkExceptionMacro( + << "CudaWarpForwardProjectionImageFilter requires a displacement field. Call SetDisplacementField()."); + } + else + { + dynamic_cast(fw.GetPointer()) + ->SetDisplacementField(reinterpret_cast( + const_cast(this->GetDisplacementField().GetPointer()))); + } +#endif + return fw; + } + + + template * = nullptr> + ForwardProjectionPointerType + InstantiateCudaWarpForwardProjection() + { + itkGenericExceptionMacro(<< "CudaWarpForwardProjectionImageFilter only available with 3D CudaImage of float."); + return nullptr; + } + + + template * = 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(bp.GetPointer()) + ->SetDisplacementField(reinterpret_cast( + const_cast(this->GetDisplacementField().GetPointer()))); + } +#endif + return bp; + } + + template * = nullptr> + BackProjectionPointerType + InstantiateCudaWarpBackProjection() + { + itkGenericExceptionMacro(<< "CudaWarpBackProjectionImageFilter only available with 3D CudaImage of float."); + return nullptr; + } + + template * = nullptr> BackProjectionPointerType InstantiateJosephBackAttenuatedProjection() diff --git a/include/rtkIterativeConeBeamReconstructionFilter.hxx b/include/rtkIterativeConeBeamReconstructionFilter.hxx index 53d22474b..66cf69384 100644 --- a/include/rtkIterativeConeBeamReconstructionFilter.hxx +++ b/include/rtkIterativeConeBeamReconstructionFilter.hxx @@ -59,6 +59,9 @@ IterativeConeBeamReconstructionFilter::Instan case (FP_ZENG): fw = InstantiateZengForwardProjection(); break; + case (FP_CUDAWARPRAYCAST): + fw = InstantiateCudaWarpForwardProjection(); + break; default: itkGenericExceptionMacro(<< "Unhandled --fp value."); } @@ -90,6 +93,9 @@ IterativeConeBeamReconstructionFilter::Instan case (BP_ZENG): bp = InstantiateZengBackProjection(); break; + case (BP_CUDAWARP): + bp = InstantiateCudaWarpBackProjection(); + break; default: itkGenericExceptionMacro(<< "Unhandled --bp value."); } diff --git a/test/rtkconjugategradientreconstructiontest.cxx b/test/rtkconjugategradientreconstructiontest.cxx index e20c9419f..937f90a4e 100644 --- a/test/rtkconjugategradientreconstructiontest.cxx +++ b/test/rtkconjugategradientreconstructiontest.cxx @@ -7,6 +7,7 @@ #ifdef USE_CUDA # include "itkCudaImage.h" #endif +#include /** * \file rtkconjugategradientreconstructiontest.cxx @@ -139,6 +140,26 @@ main(int, char **) CheckImageQuality(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; + using DVFImageType = itk::CudaImage; + auto dvfSource = rtk::ConstantImageSource::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(dvfSource->GetOutput())); + conjugategradient->SetForwardProjectionFilter(ConjugateGradientType::FP_CUDAWARPRAYCAST); + conjugategradient->SetBackProjectionFilter(ConjugateGradientType::BP_CUDAWARP); + TRY_AND_EXIT_ON_ITK_EXCEPTION(conjugategradient->Update()); + + CheckImageQuality(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;