19 #ifndef rtkFFTProjectionsConvolutionImageFilter_h 20 #define rtkFFTProjectionsConvolutionImageFilter_h 22 #include <itkImageToImageFilter.h> 23 #include <itkConceptChecking.h> 25 #include "rtkConfiguration.h" 44 template <
class TInputImage,
class TOutputImage,
class TFFTPrecision>
46 :
public itk::ImageToImageFilter<TInputImage, TOutputImage>
53 using Superclass = itk::ImageToImageFilter<TInputImage, TOutputImage>;
61 using IndexType =
typename InputImageType::IndexType;
62 using SizeType =
typename InputImageType::SizeType;
66 using FFTOutputImageType =
typename itk::Image<std::complex<TFFTPrecision>, TInputImage::ImageDimension>;
71 static constexpr
unsigned int ImageDimension = TOutputImage::ImageDimension;
89 itkGetConstMacro(GreatestPrimeFactor,
int);
96 itkGetConstMacro(TruncationCorrection,
double);
107 if (m_ZeroPadFactors != _arg)
109 m_ZeroPadFactors = _arg;
110 m_ZeroPadFactors[0] = std::max(m_ZeroPadFactors[0], 1);
111 m_ZeroPadFactors[1] = std::max(m_ZeroPadFactors[1], 1);
112 m_ZeroPadFactors[0] = std::min(m_ZeroPadFactors[0], 2);
113 m_ZeroPadFactors[1] = std::min(m_ZeroPadFactors[1], 2);
119 itkConceptMacro(ImageDimensionCheck, (itk::Concept::SameDimensionOrMinusOne<Self::InputImageDimension, 3>));
126 GenerateInputRequestedRegion()
override;
129 BeforeThreadedGenerateData()
override;
132 AfterThreadedGenerateData()
override;
135 ThreadedGenerateData(
const RegionType & outputRegionForThread, ThreadIdType threadId)
override;
141 virtual FFTInputImagePointer
142 PadInputImageRegion(
const RegionType & inputRegion);
144 GetPaddedImageRegion(
const RegionType & inputRegion);
148 PrintSelf(std::ostream & os, itk::Indent indent)
const override;
151 IsPrime(
int n)
const;
154 GreatestPrimeFactor(
int n)
const;
159 UpdateFFTProjectionsConvolutionKernel(
const SizeType size) = 0;
167 UpdateTruncationMirrorWeights();
172 int m_KernelDimension{ 1 };
183 double m_TruncationCorrection{ 0. };
185 GetTruncationCorrectionExtent();
196 int m_GreatestPrimeFactor{ 2 };
197 int m_BackupNumberOfThreads{ 1 };
202 #ifndef ITK_MANUAL_INSTANTIATION 203 # include "rtkFFTProjectionsConvolutionImageFilter.hxx" itk::SmartPointer< Self > Pointer
itk::SmartPointer< const Self > ConstPointer
Base class for 1D or 2D FFT based convolution of projections.
typename InputImageType::SizeType SizeType
itk::ImageToImageFilter< TInputImage, TOutputImage > Superclass
TInputImage InputImageType
typename InputImageType::RegionType RegionType
std::vector< TFFTPrecision > m_TruncationMirrorWeights
#define itkSetMacro(name, type)
typename FFTInputImageType::Pointer FFTInputImagePointer
typename InputImageType::IndexType IndexType
typename itk::Image< TFFTPrecision, TInputImage::ImageDimension > FFTInputImageType
FFTOutputImagePointer m_KernelFFT
virtual void SetZeroPadFactors(ZeroPadFactorsType _arg)
itk::Vector< int, 2 > ZeroPadFactorsType
typename FFTOutputImageType::Pointer FFTOutputImagePointer
ZeroPadFactorsType m_ZeroPadFactors
TOutputImage OutputImageType
typename itk::Image< std::complex< TFFTPrecision >, TInputImage::ImageDimension > FFTOutputImageType