RTK  2.6.0
Reconstruction Toolkit
rtkSimplexSpectralProjectionsDecompositionImageFilter.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 rtkSimplexSpectralProjectionsDecompositionImageFilter_h
20 #define rtkSimplexSpectralProjectionsDecompositionImageFilter_h
21 
22 #include <itkImageToImageFilter.h>
23 #include <itkAmoebaOptimizer.h>
27 
29 
30 namespace rtk
31 {
43 template <typename DecomposedProjectionsType,
44  typename MeasuredProjectionsType,
45  typename IncidentSpectrumImageType = itk::Image<float, 3>,
46  typename DetectorResponseImageType = itk::Image<float, 2>,
47  typename MaterialAttenuationsImageType = itk::Image<float, 2>>
49  : public itk::ImageToImageFilter<DecomposedProjectionsType, DecomposedProjectionsType>
50 {
51 public:
52  ITK_DISALLOW_COPY_AND_MOVE(SimplexSpectralProjectionsDecompositionImageFilter);
53 
59 
61  using InputImageType = DecomposedProjectionsType;
62  using OutputImageType = DecomposedProjectionsType;
63 
67  using DetectorResponseType = vnl_matrix<double>;
68  using MaterialAttenuationsType = vnl_matrix<double>;
70 
71 #ifndef ITK_FUTURE_LEGACY_REMOVE
72 
73  using VectorSpectrumImageType = itk::VectorImage<float, 2>;
76 #endif
77 
79  itkNewMacro(Self);
80 
82  itkOverrideGetNameOfClassMacro(SimplexSpectralProjectionsDecompositionImageFilter);
83 
85  void
86  SetInputDecomposedProjections(const DecomposedProjectionsType * DecomposedProjections);
87  typename DecomposedProjectionsType::ConstPointer
88  GetInputDecomposedProjections();
90 
92  void
93  SetInputMeasuredProjections(const MeasuredProjectionsType * SpectralProjections);
94  typename MeasuredProjectionsType::ConstPointer
95  GetInputMeasuredProjections();
97 
99  void
100  SetDetectorResponse(const DetectorResponseImageType * DetectorResponse);
101  typename DetectorResponseImageType::ConstPointer
102  GetDetectorResponse();
104 
106  void
107  SetMaterialAttenuations(const MaterialAttenuationsImageType * MaterialAttenuations);
108  typename MaterialAttenuationsImageType::ConstPointer
109  GetMaterialAttenuations();
111 
113  void
114  SetInputIncidentSpectrum(const IncidentSpectrumImageType * IncidentSpectrum);
115  void
116  SetInputSecondIncidentSpectrum(const IncidentSpectrumImageType * SecondIncidentSpectrum);
117 #ifndef ITK_FUTURE_LEGACY_REMOVE
118  void
119  SetInputIncidentSpectrum(const VectorSpectrumImageType * IncidentSpectrum);
120  void
121  SetInputSecondIncidentSpectrum(const VectorSpectrumImageType * SecondIncidentSpectrum);
122 #endif
123  typename IncidentSpectrumImageType::ConstPointer
124  GetInputIncidentSpectrum();
125  typename IncidentSpectrumImageType::ConstPointer
126  GetInputSecondIncidentSpectrum();
128 
130  itkGetMacro(NumberOfIterations, unsigned int);
131  itkSetMacro(NumberOfIterations, unsigned int);
133 
134  itkSetMacro(NumberOfEnergies, unsigned int);
135  itkGetMacro(NumberOfEnergies, unsigned int);
136 
137  itkSetMacro(NumberOfMaterials, unsigned int);
138  itkGetMacro(NumberOfMaterials, unsigned int);
139 
140  itkSetMacro(OptimizeWithRestarts, bool);
141  itkGetMacro(OptimizeWithRestarts, bool);
142 
143  itkSetMacro(Thresholds, ThresholdsType);
144  itkGetMacro(Thresholds, ThresholdsType);
145 
146  itkSetMacro(NumberOfSpectralBins, unsigned int);
147  itkGetMacro(NumberOfSpectralBins, unsigned int);
148 
149  itkSetMacro(OutputInverseCramerRaoLowerBound, bool);
150  itkGetMacro(OutputInverseCramerRaoLowerBound, bool);
151 
152  itkSetMacro(OutputFischerMatrix, bool);
153  itkGetMacro(OutputFischerMatrix, bool);
154 
155  itkSetMacro(LogTransformEachBin, bool);
156  itkGetMacro(LogTransformEachBin, bool);
157 
158  itkSetMacro(GuessInitialization, bool);
159  itkGetMacro(GuessInitialization, bool);
160 
161  itkSetMacro(IsSpectralCT, bool);
162  itkGetMacro(IsSpectralCT, bool);
163 
164 protected:
167 
168  void
169  GenerateOutputInformation() override;
170 
171  void
172  GenerateInputRequestedRegion() override;
173 
174  void
175  BeforeThreadedGenerateData() override;
176  void
177  DynamicThreadedGenerateData(const typename DecomposedProjectionsType::RegionType & outputRegionForThread) override;
178 
181  using Superclass::MakeOutput;
183  MakeOutput(DataObjectPointerArraySizeType idx) override;
184 
187  void
188  VerifyInputInformation() const override
189  {}
190 
200  bool m_IsSpectralCT; // If not, it is dual energy CT
202  unsigned int m_NumberOfIterations;
203  unsigned int m_NumberOfMaterials;
204  unsigned int m_NumberOfEnergies;
206 
207 #ifndef ITK_FUTURE_LEGACY_REMOVE
208 
209  typename FlattenVectorFilterType::Pointer m_FlattenFilter;
210  typename FlattenVectorFilterType::Pointer m_FlattenSecondFilter;
211  typename PermuteFilterType::Pointer m_PermuteFilter;
212  typename PermuteFilterType::Pointer m_PermuteSecondFilter;
213 #endif
214 
215 }; // end of class
216 
217 } // end namespace rtk
218 
219 
220 #ifndef ITK_MANUAL_INSTANTIATION
221 # include "rtkSimplexSpectralProjectionsDecompositionImageFilter.hxx"
222 #endif
223 
224 #endif
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
Decomposition of spectral projection images into material projections.
#define itkSetMacro(name, type)
Re-writes a vector image as an image.