RTK  2.6.0
Reconstruction Toolkit
rtkI0EstimationProjectionFilter.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 rtkI0EstimationProjectionFilter_h
20 #define rtkI0EstimationProjectionFilter_h
21 
22 #include <itkInPlaceImageFilter.h>
24 #include "rtkConfiguration.h"
25 
26 #include <mutex>
27 #include <vector>
28 #include <string>
29 
30 namespace rtk
31 {
43 template <class TInputImage = itk::Image<unsigned short, 3>,
44  class TOutputImage = TInputImage,
45  unsigned char bitShift = 2>
46 class ITK_TEMPLATE_EXPORT I0EstimationProjectionFilter : public itk::InPlaceImageFilter<TInputImage, TOutputImage>
47 {
48 public:
49  ITK_DISALLOW_COPY_AND_MOVE(I0EstimationProjectionFilter);
50 
56 
58  itkNewMacro(Self);
59 
61  itkOverrideGetNameOfClassMacro(I0EstimationProjectionFilter);
62 
64  using InputImageType = TInputImage;
65  using ImagePointer = typename InputImageType::Pointer;
66  using ImageConstPointer = typename InputImageType::ConstPointer;
67  using OutputImageRegionType = typename Superclass::OutputImageRegionType;
68  using InputImagePixelType = typename InputImageType::PixelType;
69 
70  itkConceptMacro(InputImagePixelTypeIsInteger, (itk::Concept::IsInteger<InputImagePixelType>));
71 
73  itkGetMacro(I0, InputImagePixelType);
74  itkGetMacro(I0fwhm, InputImagePixelType);
75  itkGetMacro(I0rls, InputImagePixelType);
77 
81  itkSetMacro(MaxPixelValue, InputImagePixelType);
82  itkGetMacro(MaxPixelValue, InputImagePixelType);
84 
86  itkSetMacro(ExpectedI0, InputImagePixelType);
87  itkGetMacro(ExpectedI0, InputImagePixelType);
89 
91  itkSetMacro(Lambda, float);
92  itkGetMacro(Lambda, float);
94 
97  itkSetMacro(Reset, bool);
98  itkGetConstMacro(Reset, bool);
99  itkBooleanMacro(Reset);
101 
104  itkSetMacro(SaveHistograms, bool);
105  itkGetConstMacro(SaveHistograms, bool);
106  itkBooleanMacro(SaveHistograms);
108 
109 protected:
111  ~I0EstimationProjectionFilter() override = default;
112 
113  void
114  BeforeThreadedGenerateData() override;
115 
116  void
117  ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread, ThreadIdType threadId) override;
118 
119  void
120  AfterThreadedGenerateData() override;
121 
122 private:
123  // Input variables
124  InputImagePixelType m_ExpectedI0; // Expected I0 value (as a result of a
125  // detector calibration)
126  InputImagePixelType m_MaxPixelValue; // Maximum encodable detector value if
127  // different from (2^16-1)
128  float m_Lambda; // RLS coefficient
129  bool m_SaveHistograms; // Save histograms in a output file
130  bool m_Reset; // Reset counters
131 
132  // Secondary inputs
133  std::vector<unsigned int>::size_type m_NBins; // Histogram size, computed
134  // from m_MaxPixelValue and bitshift
135 
136  // Main variables
137  std::vector<unsigned int> m_Histogram; // compressed (bitshifted) histogram
138  InputImagePixelType m_I0; // I0 estimate with no a priori for
139  // each new image
140  InputImagePixelType m_I0rls; // Updated RLS estimate
141  InputImagePixelType m_I0fwhm; // FWHM of the I0 mode
142 
143  // Secondary variables
144  unsigned int m_Np; // Number of previously analyzed
145  // images
146  InputImagePixelType m_Imin, m_Imax; // Define the range of consistent
147  // pixels in histogram
148  unsigned int m_DynThreshold; // Detector values with a frequency of
149  // less than dynThreshold outside
150  // min/max are discarded
151  unsigned int m_LowBound, m_HighBound; // Lower/Upper bounds of the I0 mode
152  // at half width
153 
154  std::mutex m_Mutex;
155  int m_Nsync;
157 };
158 } // end namespace rtk
159 
160 #ifndef ITK_MANUAL_INSTANTIATION
161 # include "rtkI0EstimationProjectionFilter.hxx"
162 #endif
163 
164 #endif
Estimate the I0 value from the projection histograms.
typename InputImageType::ConstPointer ImageConstPointer
TInputImage InputImageType
typename InputImageType::PixelType InputImagePixelType
#define itkSetMacro(name, type)
typename OutputImageType::RegionType OutputImageRegionType
std::vector< unsigned int >::size_type m_NBins
typename InputImageType::Pointer ImagePointer
unsigned int ThreadIdType
#define itkConceptMacro(name, concept)