RTK  2.7.0
Reconstruction Toolkit
rtkSpectralForwardModelImageFilter.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 rtkSpectralForwardModelImageFilter_h
20 #define rtkSpectralForwardModelImageFilter_h
21 
24 
26 #include <itkInPlaceImageFilter.h>
27 #include <itkCastImageFilter.h>
28 
29 namespace rtk
30 {
44 template <typename DecomposedProjectionsType,
45  typename MeasuredProjectionsType,
46  typename IncidentSpectrumImageType = itk::Image<float, 3>,
47  typename DetectorResponseImageType = itk::Image<float, 2>,
48  typename MaterialAttenuationsImageType = itk::Image<float, 2>>
49 class ITK_TEMPLATE_EXPORT SpectralForwardModelImageFilter
50  : public itk::InPlaceImageFilter<MeasuredProjectionsType, MeasuredProjectionsType>
51 {
52 public:
53  ITK_DISALLOW_COPY_AND_MOVE(SpectralForwardModelImageFilter);
54 
60 
62  using InputImageType = MeasuredProjectionsType;
63  using OutputImageType = MeasuredProjectionsType;
64 
67  using DetectorResponseType = vnl_matrix<double>;
68  using MaterialAttenuationsType = vnl_matrix<double>;
69  using DecomposedProjectionsDataType = typename DecomposedProjectionsType::PixelType::ValueType;
70  using MeasuredProjectionsDataType = typename MeasuredProjectionsType::PixelType::ValueType;
71 
72 #ifndef ITK_FUTURE_LEGACY_REMOVE
73 
74  using VectorSpectrumImageType = itk::VectorImage<float, 2>;
77 #endif
78 
80  itkNewMacro(Self);
81 
83  itkOverrideGetNameOfClassMacro(SpectralForwardModelImageFilter);
84 
86  void
87  SetInputIncidentSpectrum(const IncidentSpectrumImageType * IncidentSpectrum);
88 #ifndef ITK_FUTURE_LEGACY_REMOVE
89  void
90  SetInputIncidentSpectrum(const VectorSpectrumImageType * IncidentSpectrum);
91 #endif
92  typename IncidentSpectrumImageType::ConstPointer
93  GetInputIncidentSpectrum();
95 
97  void
98  SetInputDecomposedProjections(
99  const typename itk::ImageBase<DecomposedProjectionsType::ImageDimension> * DecomposedProjections);
100  template <unsigned int VNumberOfMaterials>
101  void
102  SetInputFixedVectorLengthDecomposedProjections(
104  DecomposedProjectionsType::ImageDimension> * DecomposedProjections);
105  typename DecomposedProjectionsType::ConstPointer
106  GetInputDecomposedProjections();
108 
110  void
111  SetInputMeasuredProjections(
112  const typename itk::ImageBase<MeasuredProjectionsType::ImageDimension> * MeasuredProjections);
113  template <unsigned int VNumberOfSpectralBins>
114  void
115  SetInputFixedVectorLengthMeasuredProjections(
117  MeasuredProjectionsType::ImageDimension> * MeasuredProjections);
118  typename MeasuredProjectionsType::ConstPointer
119  GetInputMeasuredProjections();
121 
123  void
124  SetDetectorResponse(const DetectorResponseImageType * DetectorResponse);
125  typename DetectorResponseImageType::ConstPointer
126  GetDetectorResponse();
128 
130  void
131  SetMaterialAttenuations(const MaterialAttenuationsImageType * MaterialAttenuations);
132  typename MaterialAttenuationsImageType::ConstPointer
133  GetMaterialAttenuations();
135 
136  typename DecomposedProjectionsType::ConstPointer
137  GetOutputCramerRaoLowerBound();
138 
139  typename MeasuredProjectionsType::ConstPointer
140  GetOutputVariances();
141 
142  itkSetMacro(Thresholds, ThresholdsType);
143  itkGetMacro(Thresholds, ThresholdsType);
144 
145  itkSetMacro(NumberOfSpectralBins, unsigned int);
146  itkGetMacro(NumberOfSpectralBins, unsigned int);
147 
148  itkSetMacro(NumberOfMaterials, unsigned int);
149  itkGetMacro(NumberOfMaterials, unsigned int);
150 
151  itkSetMacro(NumberOfEnergies, unsigned int);
152  itkGetMacro(NumberOfEnergies, unsigned int);
153 
154  itkSetMacro(ComputeVariances, bool);
155  itkGetMacro(ComputeVariances, bool);
156 
157  itkSetMacro(ComputeCramerRaoLowerBound, bool);
158  itkGetMacro(ComputeCramerRaoLowerBound, bool);
159 
160 protected:
162  ~SpectralForwardModelImageFilter() override = default;
163 
166  using Superclass::MakeOutput;
168  MakeOutput(DataObjectPointerArraySizeType idx) override;
169 
170  void
171  GenerateOutputInformation() override;
172 
173  void
174  GenerateInputRequestedRegion() override;
175 
176  void
177  BeforeThreadedGenerateData() override;
178  void
179  DynamicThreadedGenerateData(const typename OutputImageType::RegionType & outputRegionForThread) override;
180 
183  void
184  VerifyInputInformation() const override
185  {}
186 
189 
192  unsigned int m_NumberOfEnergies;
193 
195  unsigned int m_NumberOfIterations;
196  unsigned int m_NumberOfMaterials;
198  bool m_ComputeVariances; // Only implemented for dual energy CT
199  bool m_ComputeCramerRaoLowerBound; // Only implemented for spectral CT
200 
207 
208 #ifndef ITK_FUTURE_LEGACY_REMOVE
209 
210  typename FlattenVectorFilterType::Pointer m_FlattenFilter;
211  typename PermuteFilterType::Pointer m_PermuteFilter;
212 #endif
213 }; // end of class
214 
215 // Function to bin a detector response matrix according to given energy thresholds
216 template <typename OutputElementType, typename DetectorResponseImageType, typename ThresholdsType>
217 vnl_matrix<OutputElementType>
218 SpectralBinDetectorResponse(const DetectorResponseImageType * drm,
219  const ThresholdsType & thresholds,
220  const unsigned int numberOfEnergies);
221 
222 } // end namespace rtk
223 
224 
225 #ifndef ITK_MANUAL_INSTANTIATION
226 # include "rtkSpectralForwardModelImageFilter.hxx"
227 #endif
228 
229 #endif
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
Forward model for the decomposition of spectral projection images into material projections.
itk::ImageSource< MeasuredProjectionsType >::Pointer m_CastMeasuredProjectionsPointer
#define itkSetMacro(name, type)
typename DecomposedProjectionsType::PixelType::ValueType DecomposedProjectionsDataType
vnl_matrix< OutputElementType > SpectralBinDetectorResponse(const DetectorResponseImageType *drm, const ThresholdsType &thresholds, const unsigned int numberOfEnergies)
Re-writes a vector image as an image.
itk::ImageSource< DecomposedProjectionsType >::Pointer m_CastDecomposedProjectionsPointer
typename MeasuredProjectionsType::PixelType::ValueType MeasuredProjectionsDataType