RTK  1.4.0
Reconstruction Toolkit
rtkSingularValueThresholdImageFilter.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  * http://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 rtkSingularValueThresholdImageFilter_h
20 #define rtkSingularValueThresholdImageFilter_h
21 
22 #include <itkInPlaceImageFilter.h>
23 #include <itkVector.h>
24 
26 
27 namespace rtk
28 {
49 template< typename TInputImage,
50  typename TRealType = float,
51  typename TOutputImage = TInputImage>
53  public itk::InPlaceImageFilter< TInputImage, TOutputImage >
54 {
55 public:
56 
62 
64  itkNewMacro(Self);
65 
68 
71  typedef typename TOutputImage::PixelType OutputPixelType;
72  typedef typename TInputImage::PixelType InputPixelType;
73 
75  typedef TInputImage InputImageType;
76  typedef TOutputImage OutputImageType;
77  typedef typename InputImageType::Pointer InputImagePointer;
78  typedef typename OutputImageType::Pointer OutputImagePointer;
79 
81  itkStaticConstMacro(ImageDimension, unsigned int,
82  TOutputImage::ImageDimension);
83 
85  itkStaticConstMacro(VectorDimension, unsigned int,
86  InputPixelType::Dimension);
87 
89  typedef TRealType RealType;
92 
94  typedef typename Superclass::OutputImageRegionType OutputImageRegionType;
95 
96 #ifdef ITK_USE_CONCEPT_CHECKING
97 
98  itkConceptMacro( InputHasNumericTraitsCheck,
100  itkConceptMacro( RealTypeHasNumericTraitsCheck,
102 
104 #endif
105 
106  itkGetMacro(Threshold, TRealType)
107  itkSetMacro(Threshold, TRealType)
108 
109 protected:
111  virtual ~SingularValueThresholdImageFilter() ITK_OVERRIDE {}
112 
113  void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
114  itk::ThreadIdType threadId) ITK_OVERRIDE;
115 
117  const itk::ImageRegionSplitterBase* GetImageRegionSplitter(void) const ITK_OVERRIDE;
119 
120 private:
121  TRealType m_Threshold;
122 
123  SingularValueThresholdImageFilter(const Self &); //purposely not implemented
124  void operator=(const Self &); //purposely not implemented
125 };
126 } // end namespace itk
127 
128 #ifndef ITK_MANUAL_INSTANTIATION
129 #include "rtkSingularValueThresholdImageFilter.hxx"
130 #endif
131 
132 #endif
itk::ImageToImageFilter< TInputImage, TOutputImage > Superclass
itk::Vector< TRealType, InputPixelType::Dimension > RealVectorType
itk::Image< RealVectorType, TInputImage::ImageDimension > RealVectorImageType
itk::ImageRegionSplitterDirection::Pointer m_Splitter
unsigned int ThreadIdType
Performs thresholding on the singular values.
#define itkConceptMacro(name, concept)
#define itkSetMacro(name, type)