RTK  2.7.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>
26 
28 
29 namespace rtk
30 {
42 template <typename DecomposedProjectionsType,
43  typename MeasuredProjectionsType,
44  typename IncidentSpectrumImageType = itk::Image<float, 3>,
45  typename DetectorResponseImageType = itk::Image<float, 2>,
46  typename MaterialAttenuationsImageType = itk::Image<float, 2>>
48  : public itk::ImageToImageFilter<DecomposedProjectionsType, DecomposedProjectionsType>
49 {
50 public:
51  ITK_DISALLOW_COPY_AND_MOVE(SimplexSpectralProjectionsDecompositionImageFilter);
52 
58 
60  using InputImageType = DecomposedProjectionsType;
61  using OutputImageType = DecomposedProjectionsType;
62 
66  using DetectorResponseType = vnl_matrix<double>;
67  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(SimplexSpectralProjectionsDecompositionImageFilter);
84 
86  void
87  SetInputDecomposedProjections(
88  const typename itk::ImageBase<DecomposedProjectionsType::ImageDimension> * DecomposedProjections);
89  template <unsigned int VNumberOfMaterials>
90  void
91  SetInputFixedVectorLengthDecomposedProjections(
93  DecomposedProjectionsType::ImageDimension> * DecomposedProjections);
94  typename DecomposedProjectionsType::ConstPointer
95  GetInputDecomposedProjections();
97 
99  void
100  SetInputMeasuredProjections(
101  const typename itk::ImageBase<MeasuredProjectionsType::ImageDimension> * MeasuredProjections);
102  template <unsigned int VNumberOfSpectralBins>
103  void
104  SetInputFixedVectorLengthMeasuredProjections(
106  MeasuredProjectionsType::ImageDimension> * MeasuredProjections);
107  typename MeasuredProjectionsType::ConstPointer
108  GetInputMeasuredProjections();
110 
111 
113  void
114  SetDetectorResponse(const DetectorResponseImageType * DetectorResponse);
115  typename DetectorResponseImageType::ConstPointer
116  GetDetectorResponse();
118 
120  void
121  SetMaterialAttenuations(const MaterialAttenuationsImageType * MaterialAttenuations);
122  typename MaterialAttenuationsImageType::ConstPointer
123  GetMaterialAttenuations();
125 
127  void
128  SetInputIncidentSpectrum(const IncidentSpectrumImageType * IncidentSpectrum);
129 #ifndef ITK_FUTURE_LEGACY_REMOVE
130  void
131  SetInputIncidentSpectrum(const VectorSpectrumImageType * IncidentSpectrum);
132 #endif
133  typename IncidentSpectrumImageType::ConstPointer
134  GetInputIncidentSpectrum();
136 
138  itkGetMacro(NumberOfIterations, unsigned int);
139  itkSetMacro(NumberOfIterations, unsigned int);
141 
142  itkSetMacro(NumberOfEnergies, unsigned int);
143  itkGetMacro(NumberOfEnergies, unsigned int);
144 
145  itkSetMacro(NumberOfMaterials, unsigned int);
146  itkGetMacro(NumberOfMaterials, unsigned int);
147 
148  itkSetMacro(OptimizeWithRestarts, bool);
149  itkGetMacro(OptimizeWithRestarts, bool);
150 
151  itkSetMacro(Thresholds, ThresholdsType);
152  itkGetMacro(Thresholds, ThresholdsType);
153 
154  itkSetMacro(NumberOfSpectralBins, unsigned int);
155  itkGetMacro(NumberOfSpectralBins, unsigned int);
156 
157  itkSetMacro(OutputInverseCramerRaoLowerBound, bool);
158  itkGetMacro(OutputInverseCramerRaoLowerBound, bool);
159 
160  itkSetMacro(OutputFischerMatrix, bool);
161  itkGetMacro(OutputFischerMatrix, bool);
162 
163  itkSetMacro(LogTransformEachBin, bool);
164  itkGetMacro(LogTransformEachBin, bool);
165 
166  itkSetMacro(GuessInitialization, bool);
167  itkGetMacro(GuessInitialization, bool);
168 
169 protected:
172 
173  void
174  GenerateOutputInformation() override;
175 
176  void
177  GenerateInputRequestedRegion() override;
178 
179  void
180  BeforeThreadedGenerateData() override;
181  void
182  DynamicThreadedGenerateData(const typename DecomposedProjectionsType::RegionType & outputRegionForThread) override;
183 
186  using Superclass::MakeOutput;
188  MakeOutput(DataObjectPointerArraySizeType idx) override;
189 
192  void
193  VerifyInputInformation() const override
194  {}
195 
206  unsigned int m_NumberOfIterations;
207  unsigned int m_NumberOfMaterials;
208  unsigned int m_NumberOfEnergies;
210 
217 
218 #ifndef ITK_FUTURE_LEGACY_REMOVE
219 
220  typename FlattenVectorFilterType::Pointer m_FlattenFilter;
221  typename PermuteFilterType::Pointer m_PermuteFilter;
222 #endif
223 }; // end of class
224 
225 } // end namespace rtk
226 
227 
228 #ifndef ITK_MANUAL_INSTANTIATION
229 # include "rtkSimplexSpectralProjectionsDecompositionImageFilter.hxx"
230 #endif
231 
232 #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.
typename DecomposedProjectionsType::PixelType::ValueType DecomposedProjectionsDataType
itk::ImageSource< DecomposedProjectionsType >::Pointer m_CastDecomposedProjectionsPointer