OTB  10.0.0
Orfeo Toolbox
otbClampImageFilter.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2024 Centre National d'Etudes Spatiales (CNES)
3  *
4  * This file is part of Orfeo Toolbox
5  *
6  * https://www.orfeo-toolbox.org/
7  *
8  * Licensed under the Apache License, Version 2.0 (the "License");
9  * you may not use this file except in compliance with the License.
10  * You may obtain a copy of the License at
11  *
12  * http://www.apache.org/licenses/LICENSE-2.0
13  *
14  * Unless required by applicable law or agreed to in writing, software
15  * distributed under the License is distributed on an "AS IS" BASIS,
16  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
17  * See the License for the specific language governing permissions and
18  * limitations under the License.
19  */
20 
21 #ifndef otbClampImageFilter_h
22 #define otbClampImageFilter_h
23 
24 #include "otbConvertTypeFunctor.h"
25 #include "itkUnaryFunctorImageFilter.h"
26 
27 namespace otb
28 {
29 
47 template <class TInputImage, class TOutputImage = TInputImage>
48 class ITK_EXPORT ClampImageFilter
49 : public itk::ImageToImageFilter<TInputImage, TOutputImage>
50 {
51 public:
54  using Superclass = itk::ImageToImageFilter<TInputImage, TOutputImage>;
56  using Pointer = itk::SmartPointer<Self>;
57  using ConstPointer = itk::SmartPointer<const Self>;
58 
60  itkNewMacro(Self);
61 
63  itkTypeMacro(ClampImageFilter, OSEF);
64 
65 
67  using InputImageType = TInputImage;
68  using InputImageRegionType = typename InputImageType::RegionType;
69  using InputImagePixelType = typename InputImageType::PixelType;
70 
72  using OutputImageType = TOutputImage;
73  using OutputImageRegionType = typename OutputImageType::RegionType;
74  using OutputImagePixelType = typename OutputImageType::PixelType;
75  using OutputInternalPixelType = typename OutputImageType::InternalPixelType;
76  using OutputPixelValueType = typename itk::NumericTraits<OutputInternalPixelType>::ValueType;
77 
79  void ClampAbove(const OutputPixelValueType& thresh);
80 
82  void ClampBelow(const OutputPixelValueType& thresh);
83 
85  void ClampOutside(const OutputPixelValueType& lower, const OutputPixelValueType& upper);
86 
88  void SetThresholds(OutputPixelValueType lowerVal, OutputPixelValueType upperVal);
89  itkGetConstMacro(Lower, OutputPixelValueType);
90  itkGetConstMacro(Upper, OutputPixelValueType);
92 
93 
94 protected:
96  ~ClampImageFilter() override = default;
97  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
98 
99  void GenerateOutputInformation(void) override
100  {
101  Superclass::GenerateOutputInformation();
102  unsigned int sizeIn = this->GetInput()->GetNumberOfComponentsPerPixel();
103  this->GetFunctor().SetInputComponents(sizeIn);
104  this->GetOutput()->SetNumberOfComponentsPerPixel(this->GetFunctor().GetOutputSize());
105  }
106 
107  void DynamicThreadedGenerateData(const OutputImageRegionType& outputRegionForThread) override;
108 
109  FunctorType & GetFunctor() noexcept { return m_Functor; }
110  FunctorType const& GetFunctor() const noexcept { return m_Functor; }
111 
112 private:
113  ClampImageFilter(const Self&) = delete;
114  void operator=(const Self&) = delete;
115 
119 };
120 
121 
122 } // end namespace otb
123 
124 #ifndef OTB_MANUAL_INSTANTIATION
125 #include "otbClampImageFilter.hxx"
126 #endif
127 
128 #endif
Clamp image values to be below, over, or between threshold values.
itk::SmartPointer< Self > Pointer
typename OutputImageType::InternalPixelType OutputInternalPixelType
FunctorType & GetFunctor() noexcept
ClampImageFilter(const Self &)=delete
~ClampImageFilter() override=default
typename OutputImageType::RegionType OutputImageRegionType
void operator=(const Self &)=delete
itk::ImageToImageFilter< TInputImage, TOutputImage > Superclass
FunctorType const & GetFunctor() const noexcept
typename itk::NumericTraits< OutputInternalPixelType >::ValueType OutputPixelValueType
OutputPixelValueType m_Upper
OutputPixelValueType m_Lower
typename InputImageType::RegionType InputImageRegionType
typename OutputImageType::PixelType OutputImagePixelType
typename InputImageType::PixelType InputImagePixelType
void GenerateOutputInformation(void) override
itk::SmartPointer< const Self > ConstPointer
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.