RTK  2.6.0
Reconstruction Toolkit
rtkBackwardDifferenceDivergenceImageFilter.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 rtkBackwardDifferenceDivergenceImageFilter_h
20 #define rtkBackwardDifferenceDivergenceImageFilter_h
21 
22 #include <itkImageToImageFilter.h>
23 #include <itkCastImageFilter.h>
24 
25 namespace rtk
26 {
39 template <typename TInputImage, typename TOutputImage = itk::Image<float, TInputImage::ImageDimension>>
40 class ITK_TEMPLATE_EXPORT BackwardDifferenceDivergenceImageFilter
41  : public itk::ImageToImageFilter<TInputImage, TOutputImage>
42 {
43 public:
44  ITK_DISALLOW_COPY_AND_MOVE(BackwardDifferenceDivergenceImageFilter);
45 
47  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
48 
50  using InputImageType = TInputImage;
51 
57 
59  itkNewMacro(Self);
60 
62  itkOverrideGetNameOfClassMacro(BackwardDifferenceDivergenceImageFilter);
63 
66  void
68  {
69  this->SetUseImageSpacing(true);
70  }
71 
74  void
76  {
77  this->SetUseImageSpacing(false);
78  }
79 
82  itkSetMacro(UseImageSpacing, bool);
83  itkGetConstMacro(UseImageSpacing, bool);
85 
88  void
89  SetDimensionsProcessed(bool * DimensionsProcessed);
90 
92  void
93  OverrideBoundaryCondition(itk::ImageBoundaryCondition<TInputImage> * boundaryCondition);
94 
96  using InputPixelType = typename InputImageType::PixelType;
97  using InputImageRegionType = typename InputImageType::RegionType;
98  using InputSizeType = typename InputImageType::SizeType;
100 
101 protected:
104 
105  void
106  GenerateInputRequestedRegion() override;
107 
108  void
109  BeforeThreadedGenerateData() override;
110 
111  void
112  DynamicThreadedGenerateData(const typename InputImageType::RegionType & outputRegionForThread) override;
113 
114  void
115  AfterThreadedGenerateData() override;
116 
117 private:
119  typename TInputImage::SpacingType m_InvSpacingCoeffs;
120 
121  // list of the dimensions along which the divergence has
122  // to be computed. The components on other dimensions
123  // are ignored for performance, but the gradient filter
124  // sets them to zero anyway
125  bool m_DimensionsProcessed[TInputImage::ImageDimension];
126 
127  // The default is ConstantBoundaryCondition, but this behavior sometimes needs to be overriden
129  // If so, do not perform boundary processing in AfterThreadedGenerateData
131 };
132 
133 } // namespace rtk
134 
135 #ifndef ITK_MANUAL_INSTANTIATION
136 # include "rtkBackwardDifferenceDivergenceImageFilter.hxx"
137 #endif
138 
139 #endif //__rtkBackwardDifferenceDivergenceImageFilter__
#define itkSetMacro(name, type)
itk::ImageBoundaryCondition< TInputImage, TInputImage > * m_BoundaryCondition
Computes the backward differences divergence (adjoint of the forward differences gradient) of the inp...