RTK  2.6.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 
25 
27 #include <itkInPlaceImageFilter.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 
70 #ifndef ITK_FUTURE_LEGACY_REMOVE
71 
72  using VectorSpectrumImageType = itk::VectorImage<float, 2>;
75 #endif
76 
78  itkNewMacro(Self);
79 
81  itkOverrideGetNameOfClassMacro(SpectralForwardModelImageFilter);
82 
84  void
85  SetInputIncidentSpectrum(const IncidentSpectrumImageType * IncidentSpectrum);
86  void
87  SetInputSecondIncidentSpectrum(const IncidentSpectrumImageType * SecondIncidentSpectrum);
88 #ifndef ITK_FUTURE_LEGACY_REMOVE
89  void
90  SetInputIncidentSpectrum(const VectorSpectrumImageType * IncidentSpectrum);
91  void
92  SetInputSecondIncidentSpectrum(const VectorSpectrumImageType * SecondIncidentSpectrum);
93 #endif
94  typename IncidentSpectrumImageType::ConstPointer
95  GetInputIncidentSpectrum();
96  typename IncidentSpectrumImageType::ConstPointer
97  GetInputSecondIncidentSpectrum();
99 
101  void
102  SetInputDecomposedProjections(const DecomposedProjectionsType * DecomposedProjections);
103  typename DecomposedProjectionsType::ConstPointer
104  GetInputDecomposedProjections();
106 
108  void
109  SetInputMeasuredProjections(const MeasuredProjectionsType * SpectralProjections);
110  typename MeasuredProjectionsType::ConstPointer
111  GetInputMeasuredProjections();
113 
115  void
116  SetDetectorResponse(const DetectorResponseImageType * DetectorResponse);
117  typename DetectorResponseImageType::ConstPointer
118  GetDetectorResponse();
120 
122  void
123  SetMaterialAttenuations(const MaterialAttenuationsImageType * MaterialAttenuations);
124  typename MaterialAttenuationsImageType::ConstPointer
125  GetMaterialAttenuations();
127 
128  typename DecomposedProjectionsType::ConstPointer
129  GetOutputCramerRaoLowerBound();
130 
131  typename MeasuredProjectionsType::ConstPointer
132  GetOutputVariances();
133 
134  itkSetMacro(Thresholds, ThresholdsType);
135  itkGetMacro(Thresholds, ThresholdsType);
136 
137  itkSetMacro(NumberOfSpectralBins, unsigned int);
138  itkGetMacro(NumberOfSpectralBins, unsigned int);
139 
140  itkSetMacro(NumberOfMaterials, unsigned int);
141  itkGetMacro(NumberOfMaterials, unsigned int);
142 
143  itkSetMacro(NumberOfEnergies, unsigned int);
144  itkGetMacro(NumberOfEnergies, unsigned int);
145 
146  itkSetMacro(IsSpectralCT, bool);
147  itkGetMacro(IsSpectralCT, bool);
148 
149  itkSetMacro(ComputeVariances, bool);
150  itkGetMacro(ComputeVariances, bool);
151 
152  itkSetMacro(ComputeCramerRaoLowerBound, bool);
153  itkGetMacro(ComputeCramerRaoLowerBound, bool);
154 
155 protected:
157  ~SpectralForwardModelImageFilter() override = default;
158 
161  using Superclass::MakeOutput;
163  MakeOutput(DataObjectPointerArraySizeType idx) override;
164 
165  void
166  GenerateOutputInformation() override;
167 
168  void
169  GenerateInputRequestedRegion() override;
170 
171  void
172  BeforeThreadedGenerateData() override;
173  void
174  DynamicThreadedGenerateData(const typename OutputImageType::RegionType & outputRegionForThread) override;
175 
178  void
179  VerifyInputInformation() const override
180  {}
181 
184 
187  unsigned int m_NumberOfEnergies;
188 
190  unsigned int m_NumberOfIterations;
191  unsigned int m_NumberOfMaterials;
193  bool m_IsSpectralCT; // If not, it is dual energy CT
194  bool m_ComputeVariances; // Only implemented for dual energy CT
195  bool m_ComputeCramerRaoLowerBound; // Only implemented for spectral CT
196 
197 #ifndef ITK_FUTURE_LEGACY_REMOVE
198 
199  typename FlattenVectorFilterType::Pointer m_FlattenFilter;
200  typename FlattenVectorFilterType::Pointer m_FlattenSecondFilter;
201  typename PermuteFilterType::Pointer m_PermuteFilter;
202  typename PermuteFilterType::Pointer m_PermuteSecondFilter;
203 #endif
204 
205 }; // end of class
206 
207 // Function to bin a detector response matrix according to given energy thresholds
208 template <typename OutputElementType, typename DetectorResponseImageType, typename ThresholdsType>
209 vnl_matrix<OutputElementType>
210 SpectralBinDetectorResponse(const DetectorResponseImageType * drm,
211  const ThresholdsType & thresholds,
212  const unsigned int numberOfEnergies);
213 
214 } // end namespace rtk
215 
216 
217 #ifndef ITK_MANUAL_INSTANTIATION
218 # include "rtkSpectralForwardModelImageFilter.hxx"
219 #endif
220 
221 #endif
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
Forward model for the decomposition of spectral projection images into material projections.
#define itkSetMacro(name, type)
vnl_matrix< OutputElementType > SpectralBinDetectorResponse(const DetectorResponseImageType *drm, const ThresholdsType &thresholds, const unsigned int numberOfEnergies)
Re-writes a vector image as an image.