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>
26 
27 namespace rtk
28 {
40 template <typename DecomposedProjectionsType,
41  typename MeasuredProjectionsType,
42  typename IncidentSpectrumImageType = itk::VectorImage<float, 2>,
43  typename DetectorResponseImageType = itk::Image<float, 2>,
44  typename MaterialAttenuationsImageType = itk::Image<float, 2>>
46  : public itk::ImageToImageFilter<DecomposedProjectionsType, DecomposedProjectionsType>
47 {
48 public:
49  ITK_DISALLOW_COPY_AND_MOVE(SimplexSpectralProjectionsDecompositionImageFilter);
50 
56 
58  using InputImageType = DecomposedProjectionsType;
59  using OutputImageType = DecomposedProjectionsType;
60 
64  using DetectorResponseType = vnl_matrix<double>;
65  using MaterialAttenuationsType = vnl_matrix<double>;
67 
69  itkNewMacro(Self);
70 
72  itkOverrideGetNameOfClassMacro(SimplexSpectralProjectionsDecompositionImageFilter);
73 
75  void
76  SetInputDecomposedProjections(const DecomposedProjectionsType * DecomposedProjections);
77  typename DecomposedProjectionsType::ConstPointer
78  GetInputDecomposedProjections();
80 
82  void
83  SetInputMeasuredProjections(const MeasuredProjectionsType * SpectralProjections);
84  typename MeasuredProjectionsType::ConstPointer
85  GetInputMeasuredProjections();
87 
89  void
90  SetDetectorResponse(const DetectorResponseImageType * DetectorResponse);
91  typename DetectorResponseImageType::ConstPointer
92  GetDetectorResponse();
94 
96  void
97  SetMaterialAttenuations(const MaterialAttenuationsImageType * MaterialAttenuations);
98  typename MaterialAttenuationsImageType::ConstPointer
99  GetMaterialAttenuations();
101 
103  void
104  SetInputIncidentSpectrum(const IncidentSpectrumImageType * IncidentSpectrum);
105  void
106  SetInputSecondIncidentSpectrum(const IncidentSpectrumImageType * SecondIncidentSpectrum);
107  typename IncidentSpectrumImageType::ConstPointer
108  GetInputIncidentSpectrum();
109  typename IncidentSpectrumImageType::ConstPointer
110  GetInputSecondIncidentSpectrum();
112 
114  itkGetMacro(NumberOfIterations, unsigned int);
115  itkSetMacro(NumberOfIterations, unsigned int);
117 
118  itkSetMacro(NumberOfEnergies, unsigned int);
119  itkGetMacro(NumberOfEnergies, unsigned int);
120 
121  itkSetMacro(NumberOfMaterials, unsigned int);
122  itkGetMacro(NumberOfMaterials, unsigned int);
123 
124  itkSetMacro(OptimizeWithRestarts, bool);
125  itkGetMacro(OptimizeWithRestarts, bool);
126 
127  itkSetMacro(Thresholds, ThresholdsType);
128  itkGetMacro(Thresholds, ThresholdsType);
129 
130  itkSetMacro(NumberOfSpectralBins, unsigned int);
131  itkGetMacro(NumberOfSpectralBins, unsigned int);
132 
133  itkSetMacro(OutputInverseCramerRaoLowerBound, bool);
134  itkGetMacro(OutputInverseCramerRaoLowerBound, bool);
135 
136  itkSetMacro(OutputFischerMatrix, bool);
137  itkGetMacro(OutputFischerMatrix, bool);
138 
139  itkSetMacro(LogTransformEachBin, bool);
140  itkGetMacro(LogTransformEachBin, bool);
141 
142  itkSetMacro(GuessInitialization, bool);
143  itkGetMacro(GuessInitialization, bool);
144 
145  itkSetMacro(IsSpectralCT, bool);
146  itkGetMacro(IsSpectralCT, bool);
147 
148 protected:
151 
152  void
153  GenerateOutputInformation() override;
154 
155  void
156  GenerateInputRequestedRegion() override;
157 
158  void
159  BeforeThreadedGenerateData() override;
160  void
161  DynamicThreadedGenerateData(const typename DecomposedProjectionsType::RegionType & outputRegionForThread) override;
162 
165  using Superclass::MakeOutput;
167  MakeOutput(DataObjectPointerArraySizeType idx) override;
168 
171  void
172  VerifyInputInformation() const override
173  {}
174 
184  bool m_IsSpectralCT; // If not, it is dual energy CT
186  unsigned int m_NumberOfIterations;
187  unsigned int m_NumberOfMaterials;
188  unsigned int m_NumberOfEnergies;
190 
191 }; // end of class
192 
193 } // end namespace rtk
194 
195 
196 #ifndef ITK_MANUAL_INSTANTIATION
197 # include "rtkSimplexSpectralProjectionsDecompositionImageFilter.hxx"
198 #endif
199 
200 #endif
Superclass::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
Decomposition of spectral projection images into material projections.
#define itkSetMacro(name, type)
TOutputImage OutputImageType