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 
24 
25 #include <itkInPlaceImageFilter.h>
26 
27 namespace rtk
28 {
42 template <typename DecomposedProjectionsType,
43  typename MeasuredProjectionsType,
44  typename IncidentSpectrumImageType = itk::VectorImage<float, 2>,
45  typename DetectorResponseImageType = itk::Image<float, 2>,
46  typename MaterialAttenuationsImageType = itk::Image<float, 2>>
47 class ITK_TEMPLATE_EXPORT SpectralForwardModelImageFilter
48  : public itk::InPlaceImageFilter<MeasuredProjectionsType, MeasuredProjectionsType>
49 {
50 public:
51  ITK_DISALLOW_COPY_AND_MOVE(SpectralForwardModelImageFilter);
52 
58 
60  using InputImageType = MeasuredProjectionsType;
61  using OutputImageType = MeasuredProjectionsType;
62 
65  using DetectorResponseType = vnl_matrix<double>;
66  using MaterialAttenuationsType = vnl_matrix<double>;
67 
69  itkNewMacro(Self);
70 
72  itkOverrideGetNameOfClassMacro(SpectralForwardModelImageFilter);
73 
75  void
76  SetInputIncidentSpectrum(const IncidentSpectrumImageType * IncidentSpectrum);
77  void
78  SetInputSecondIncidentSpectrum(const IncidentSpectrumImageType * SecondIncidentSpectrum);
79  typename IncidentSpectrumImageType::ConstPointer
80  GetInputIncidentSpectrum();
81  typename IncidentSpectrumImageType::ConstPointer
82  GetInputSecondIncidentSpectrum();
84 
86  void
87  SetInputDecomposedProjections(const DecomposedProjectionsType * DecomposedProjections);
88  typename DecomposedProjectionsType::ConstPointer
89  GetInputDecomposedProjections();
91 
93  void
94  SetInputMeasuredProjections(const MeasuredProjectionsType * SpectralProjections);
95  typename MeasuredProjectionsType::ConstPointer
96  GetInputMeasuredProjections();
98 
100  void
101  SetDetectorResponse(const DetectorResponseImageType * DetectorResponse);
102  typename DetectorResponseImageType::ConstPointer
103  GetDetectorResponse();
105 
107  void
108  SetMaterialAttenuations(const MaterialAttenuationsImageType * MaterialAttenuations);
109  typename MaterialAttenuationsImageType::ConstPointer
110  GetMaterialAttenuations();
112 
113  typename DecomposedProjectionsType::ConstPointer
114  GetOutputCramerRaoLowerBound();
115 
116  typename MeasuredProjectionsType::ConstPointer
117  GetOutputVariances();
118 
119  itkSetMacro(Thresholds, ThresholdsType);
120  itkGetMacro(Thresholds, ThresholdsType);
121 
122  itkSetMacro(NumberOfSpectralBins, unsigned int);
123  itkGetMacro(NumberOfSpectralBins, unsigned int);
124 
125  itkSetMacro(NumberOfMaterials, unsigned int);
126  itkGetMacro(NumberOfMaterials, unsigned int);
127 
128  itkSetMacro(NumberOfEnergies, unsigned int);
129  itkGetMacro(NumberOfEnergies, unsigned int);
130 
131  itkSetMacro(IsSpectralCT, bool);
132  itkGetMacro(IsSpectralCT, bool);
133 
134  itkSetMacro(ComputeVariances, bool);
135  itkGetMacro(ComputeVariances, bool);
136 
137  itkSetMacro(ComputeCramerRaoLowerBound, bool);
138  itkGetMacro(ComputeCramerRaoLowerBound, bool);
139 
140 protected:
142  ~SpectralForwardModelImageFilter() override = default;
143 
146  using Superclass::MakeOutput;
148  MakeOutput(DataObjectPointerArraySizeType idx) override;
149 
150  void
151  GenerateOutputInformation() override;
152 
153  void
154  GenerateInputRequestedRegion() override;
155 
156  void
157  BeforeThreadedGenerateData() override;
158  void
159  DynamicThreadedGenerateData(const typename OutputImageType::RegionType & outputRegionForThread) override;
160 
163  void
164  VerifyInputInformation() const override
165  {}
166 
169 
172  unsigned int m_NumberOfEnergies;
173 
175  unsigned int m_NumberOfIterations;
176  unsigned int m_NumberOfMaterials;
178  bool m_IsSpectralCT; // If not, it is dual energy CT
179  bool m_ComputeVariances; // Only implemented for dual energy CT
180  bool m_ComputeCramerRaoLowerBound; // Only implemented for spectral CT
181 
182 }; // end of class
183 
184 // Function to bin a detector response matrix according to given energy thresholds
185 template <typename OutputElementType, typename DetectorResponseImageType, typename ThresholdsType>
186 vnl_matrix<OutputElementType>
187 SpectralBinDetectorResponse(const DetectorResponseImageType * drm,
188  const ThresholdsType & thresholds,
189  const unsigned int numberOfEnergies);
190 
191 } // end namespace rtk
192 
193 
194 #ifndef ITK_MANUAL_INSTANTIATION
195 # include "rtkSpectralForwardModelImageFilter.hxx"
196 #endif
197 
198 #endif
Superclass::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
Forward model for the decomposition of spectral projection images into material projections.
TInputImage InputImageType
#define itkSetMacro(name, type)
TOutputImage OutputImageType
vnl_matrix< OutputElementType > SpectralBinDetectorResponse(const DetectorResponseImageType *drm, const ThresholdsType &thresholds, const unsigned int numberOfEnergies)