RTK  2.6.0
Reconstruction Toolkit
rtkJosephForwardProjectionImageFilter.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 rtkJosephForwardProjectionImageFilter_h
20 #define rtkJosephForwardProjectionImageFilter_h
21 
22 #include "rtkConfiguration.h"
24 #include "rtkMacro.h"
25 #include <itkPixelTraits.h>
26 
28 
29 #include <itkVectorImage.h>
30 namespace rtk
31 {
32 namespace Functor
33 {
42 template <class TInput, class TCoordinateType, class TOutput = TInput>
43 class ITK_TEMPLATE_EXPORT InterpolationWeightMultiplication
44 {
45 public:
48  bool
50  {
51  return false;
52  }
53  bool
55  {
56  return !(*this != other);
57  }
58 
59  inline TOutput
60  operator()(const ThreadIdType itkNotUsed(threadId),
61  const double itkNotUsed(stepLengthInVoxel),
62  const TCoordinateType weight,
63  const TInput * p,
64  const int i) const
65  {
66  return weight * p[i];
67  }
68 };
69 
77 template <class TInput, class TOutput>
78 class ITK_TEMPLATE_EXPORT SumAlongRay
79 {
80 public:
82 
83  SumAlongRay() = default;
84  ~SumAlongRay() = default;
85  bool
86  operator!=(const SumAlongRay &) const
87  {
88  return false;
89  }
90  bool
91  operator==(const SumAlongRay & other) const
92  {
93  return !(*this != other);
94  }
95 
96  inline void
97  operator()(const ThreadIdType itkNotUsed(threadId),
98  TOutput & sumValue,
99  const TInput volumeValue,
100  const VectorType & itkNotUsed(stepInMM))
101  {
102  sumValue += static_cast<TOutput>(volumeValue);
103  }
104 };
105 
113 template <class TInput, class TOutput>
114 class ITK_TEMPLATE_EXPORT ProjectedValueAccumulation
115 {
116 public:
118 
119  ProjectedValueAccumulation() = default;
120  ~ProjectedValueAccumulation() = default;
121  bool
123  {
124  return false;
125  }
126  bool
128  {
129  return !(*this != other);
130  }
131 
132  inline void
133  operator()(const ThreadIdType itkNotUsed(threadId),
134  const TInput & input,
135  TOutput & output,
136  const TOutput & rayCastValue,
137  const VectorType & stepInMM,
138  const VectorType & itkNotUsed(source),
139  const VectorType & itkNotUsed(sourceToPixel),
140  const VectorType & itkNotUsed(nearestPoint),
141  const VectorType & itkNotUsed(farthestPoint)) const
142  {
143  output = input + rayCastValue * stepInMM.GetNorm();
144  }
145 };
146 
147 } // end namespace Functor
148 
149 
165 template <class TInputImage,
166  class TOutputImage,
167  class TInterpolationWeightMultiplication = Functor::InterpolationWeightMultiplication<
168  typename TInputImage::PixelType,
170  class TProjectedValueAccumulation =
171  Functor::ProjectedValueAccumulation<typename TInputImage::PixelType, typename TOutputImage::PixelType>,
172  class TSumAlongRay = Functor::SumAlongRay<typename TInputImage::PixelType, typename TOutputImage::PixelType>>
173 class ITK_TEMPLATE_EXPORT JosephForwardProjectionImageFilter
174  : public ForwardProjectionImageFilter<TInputImage, TOutputImage>
175 {
176 public:
177  ITK_DISALLOW_COPY_AND_MOVE(JosephForwardProjectionImageFilter);
178 
184  using InputPixelType = typename TInputImage::PixelType;
185  using OutputPixelType = typename TOutputImage::PixelType;
186  using OutputImageRegionType = typename TOutputImage::RegionType;
187  using CoordinateType = double;
191 
193  itkNewMacro(Self);
194 
196  itkOverrideGetNameOfClassMacro(JosephForwardProjectionImageFilter);
197 
199  TInterpolationWeightMultiplication &
201  {
202  return m_InterpolationWeightMultiplication;
203  }
204  const TInterpolationWeightMultiplication &
206  {
207  return m_InterpolationWeightMultiplication;
208  }
209  void
210  SetInterpolationWeightMultiplication(const TInterpolationWeightMultiplication & _arg)
211  {
212  if (m_InterpolationWeightMultiplication != _arg)
213  {
214  m_InterpolationWeightMultiplication = _arg;
215  this->Modified();
216  }
217  }
219 
222  TProjectedValueAccumulation &
224  {
225  return m_ProjectedValueAccumulation;
226  }
227  const TProjectedValueAccumulation &
229  {
230  return m_ProjectedValueAccumulation;
231  }
232  void
233  SetProjectedValueAccumulation(const TProjectedValueAccumulation & _arg)
234  {
235  if (m_ProjectedValueAccumulation != _arg)
236  {
237  m_ProjectedValueAccumulation = _arg;
238  this->Modified();
239  }
240  }
242 
244  TSumAlongRay &
246  {
247  return m_SumAlongRay;
248  }
249  const TSumAlongRay &
251  {
252  return m_SumAlongRay;
253  }
254  void
255  SetSumAlongRay(const TSumAlongRay & _arg)
256  {
257  if (m_SumAlongRay != _arg)
258  {
259  m_SumAlongRay = _arg;
260  this->Modified();
261  }
262  }
264 
268  void
269  SetInferiorClipImage(const TClipImageType * inferiorClipImage)
270  {
271  // Process object is not const-correct so the const casting is required.
272  this->SetInput("InferiorClipImage", const_cast<TClipImageType *>(inferiorClipImage));
273  }
274  typename TClipImageType::ConstPointer
276  {
277  return static_cast<const TClipImageType *>(this->itk::ProcessObject::GetInput("InferiorClipImage"));
278  }
280 
284  void
285  SetSuperiorClipImage(const TClipImageType * superiorClipImage)
286  {
287  // Process object is not const-correct so the const casting is required.
288  this->SetInput("SuperiorClipImage", const_cast<TClipImageType *>(superiorClipImage));
289  }
290  typename TClipImageType::ConstPointer
292  {
293  return static_cast<const TClipImageType *>(this->itk::ProcessObject::GetInput("SuperiorClipImage"));
294  }
296 
300  itkGetMacro(InferiorClip, double);
301  itkSetMacro(InferiorClip, double);
302  itkGetMacro(SuperiorClip, double);
303  itkSetMacro(SuperiorClip, double);
305 
306 protected:
308  ~JosephForwardProjectionImageFilter() override = default;
309 
311  void
312  GenerateInputRequestedRegion() override;
313 
314  void
315  ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread, ThreadIdType threadId) override;
316 
319  void
320  VerifyInputInformation() const override;
321 
322  inline OutputPixelType
323  BilinearInterpolation(const ThreadIdType threadId,
324  const double stepLengthInVoxel,
325  const InputPixelType * pxiyi,
326  const InputPixelType * pxsyi,
327  const InputPixelType * pxiys,
328  const InputPixelType * pxsys,
329  const double x,
330  const double y,
331  const int ox,
332  const int oy);
333 
334  inline OutputPixelType
335  BilinearInterpolationOnBorders(const ThreadIdType threadId,
336  const double stepLengthInVoxel,
337  const InputPixelType * pxiyi,
338  const InputPixelType * pxsyi,
339  const InputPixelType * pxiys,
340  const InputPixelType * pxsys,
341  const double x,
342  const double y,
343  const int ox,
344  const int oy,
345  const double minx,
346  const double miny,
347  const double maxx,
348  const double maxy);
349 
350 private:
351  // Functors
352  TInterpolationWeightMultiplication m_InterpolationWeightMultiplication;
353  TProjectedValueAccumulation m_ProjectedValueAccumulation;
354  TSumAlongRay m_SumAlongRay;
355  double m_InferiorClip{ 0. };
356  double m_SuperiorClip{ 1. };
357 };
358 
359 } // end namespace rtk
360 
361 #ifndef ITK_MANUAL_INSTANTIATION
362 # include "rtkJosephForwardProjectionImageFilter.hxx"
363 #endif
364 
365 #endif
void SetInterpolationWeightMultiplication(const TInterpolationWeightMultiplication &_arg)
Base class for forward projection, i.e. accumulation along x-ray lines.
TOutput operator()(const ThreadIdType, const double, const TCoordinateType weight, const TInput *p, const int i) const
void operator()(const ThreadIdType, TOutput &sumValue, const TInput volumeValue, const VectorType &)
Function to accumulate the ray casting on the projection.
const TProjectedValueAccumulation & GetProjectedValueAccumulation() const
Function to compute the attenuation correction on the projection.
bool operator==(const SumAlongRay &other) const
DataObject * GetInput(const DataObjectIdentifierType &key)
const TInterpolationWeightMultiplication & GetInterpolationWeightMultiplication() const
bool operator!=(const SumAlongRay &) const
typename TPixelType::ValueType ValueType
#define itkSetMacro(name, type)
bool operator!=(const ProjectedValueAccumulation &) const
Function to multiply the interpolation weights with the projected volume values.
typename OutputImageType::RegionType OutputImageRegionType
bool operator==(const InterpolationWeightMultiplication &other) const
void SetProjectedValueAccumulation(const TProjectedValueAccumulation &_arg)
TInterpolationWeightMultiplication m_InterpolationWeightMultiplication
bool operator!=(const InterpolationWeightMultiplication &) const
void SetSuperiorClipImage(const TClipImageType *superiorClipImage)
bool operator==(const ProjectedValueAccumulation &other) const
TInterpolationWeightMultiplication & GetInterpolationWeightMultiplication()
unsigned int ThreadIdType
void operator()(const ThreadIdType, const TInput &input, TOutput &output, const TOutput &rayCastValue, const VectorType &stepInMM, const VectorType &, const VectorType &, const VectorType &, const VectorType &) const
void SetInferiorClipImage(const TClipImageType *inferiorClipImage)