RTK  2.7.0
Reconstruction Toolkit
rtkDePierroRegularizationImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright RTK Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * https://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 
19 #ifndef rtkDePierroRegularizationImageFilter_h
20 #define rtkDePierroRegularizationImageFilter_h
21 
22 #include <itkMultiplyImageFilter.h>
23 #include <itkSubtractImageFilter.h>
24 #include <itkDivideImageFilter.h>
25 #include <itkImageKernelOperator.h>
26 #include <itkNeighborhoodOperatorImageFilter.h>
27 #include <itkConstantBoundaryCondition.h>
28 #include "rtkConstantImageSource.h"
29 
30 namespace rtk
31 {
32 
74 template <class TInputImage, class TOutputImage = TInputImage>
75 class ITK_TEMPLATE_EXPORT DePierroRegularizationImageFilter : public itk::ImageToImageFilter<TInputImage, TOutputImage>
76 {
77 public:
78  ITK_DISALLOW_COPY_AND_MOVE(DePierroRegularizationImageFilter);
79 
82  using Superclass = itk::ImageToImageFilter<TOutputImage, TOutputImage>;
83  using Pointer = itk::SmartPointer<Self>;
84  using ConstPointer = itk::SmartPointer<const Self>;
85 
87  using InputImageType = TInputImage;
88  using InputImagePointerType = typename TInputImage::Pointer;
89  using OutputImageType = TOutputImage;
90  using InputPixelType = typename TInputImage::PixelType;
91 
93  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
94 
96  using MultiplyImageFilterType = itk::MultiplyImageFilter<InputImageType, InputImageType>;
97  using MultpiplyImageFilterPointerType = typename MultiplyImageFilterType::Pointer;
100  using SubtractImageFilterType = itk::SubtractImageFilter<InputImageType, InputImageType>;
101  using SubtractImageFilterPointerType = typename SubtractImageFilterType::Pointer;
102  using ImageKernelOperatorType = itk::ImageKernelOperator<InputPixelType, InputImageDimension>;
103  using NOIFType = itk::NeighborhoodOperatorImageFilter<InputImageType, InputImageType>;
104  using NOIFPointerType = typename NOIFType::Pointer;
105  using CustomBinaryFilterType = itk::BinaryGeneratorImageFilter<InputImageType, InputImageType, OutputImageType>;
106  using CustomBinaryFilterPointerType = typename CustomBinaryFilterType::Pointer;
107 
109  using BoundaryCondition = itk::ConstantBoundaryCondition<InputImageType>;
110 
112  itkNewMacro(Self);
113 
115  itkOverrideGetNameOfClassMacro(DePierroRegularizationImageFilter);
116 
118  itkGetMacro(Beta, double);
119  itkSetMacro(Beta, double);
121 
122 protected:
124  ~DePierroRegularizationImageFilter() override = default;
125 
126  void
127  GenerateInputRequestedRegion() override;
128 
129  void
130  GenerateOutputInformation() override;
131 
132  void
133  GenerateData() override;
134 
144 
145 private:
146  double m_Beta{ 0.01 };
147 
148 }; // end of class
149 
150 } // end namespace rtk
151 
152 #ifndef ITK_MANUAL_INSTANTIATION
153 # include "rtkDePierroRegularizationImageFilter.hxx"
154 #endif
155 
156 #endif
itk::ConstantBoundaryCondition< InputImageType > BoundaryCondition
Generate an n-dimensional image with constant pixel values.
#define itkSetMacro(name, type)
typename ConstantVolumeSourceType::Pointer ConstantVolumeSourcePointerType
typename MultiplyImageFilterType::Pointer MultpiplyImageFilterPointerType
Implements a regularization for MLEM/OSEM reconstruction.
itk::ImageToImageFilter< TOutputImage, TOutputImage > Superclass
itk::SmartPointer< Self > Pointer
itk::NeighborhoodOperatorImageFilter< InputImageType, InputImageType > NOIFType
typename CustomBinaryFilterType::Pointer CustomBinaryFilterPointerType
itk::BinaryGeneratorImageFilter< InputImageType, InputImageType, OutputImageType > CustomBinaryFilterType
typename SubtractImageFilterType::Pointer SubtractImageFilterPointerType
itk::ImageKernelOperator< InputPixelType, InputImageDimension > ImageKernelOperatorType
itk::MultiplyImageFilter< InputImageType, InputImageType > MultiplyImageFilterType
itk::SubtractImageFilter< InputImageType, InputImageType > SubtractImageFilterType