ITK  5.4.0
Insight Toolkit
itkBinaryThresholdSpatialFunction.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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 #ifndef itkBinaryThresholdSpatialFunction_h
19 #define itkBinaryThresholdSpatialFunction_h
20 
21 #include "itkSpatialFunction.h"
22 #include "itkImageBase.h"
23 
24 namespace itk
25 {
41 template <typename TFunction>
42 class ITK_TEMPLATE_EXPORT BinaryThresholdSpatialFunction
43  : public SpatialFunction<bool, TFunction::ImageDimension, typename TFunction::InputType>
44 {
45 public:
46  ITK_DISALLOW_COPY_AND_MOVE(BinaryThresholdSpatialFunction);
47 
51 
54 
56  itkOverrideGetNameOfClassMacro(BinaryThresholdSpatialFunction);
57 
59  itkNewMacro(Self);
60 
62  using typename Superclass::OutputType;
63 
65  using InputType = typename TFunction::InputType;
66 
68  using FunctionType = TFunction;
69 
71  using FunctionOutputType = typename TFunction::OutputType;
72 
74  itkSetMacro(LowerThreshold, FunctionOutputType);
75  itkGetConstReferenceMacro(LowerThreshold, FunctionOutputType);
79  itkSetMacro(UpperThreshold, FunctionOutputType);
80  itkGetConstReferenceMacro(UpperThreshold, FunctionOutputType);
84  itkSetObjectMacro(Function, FunctionType);
85  itkGetModifiableObjectMacro(Function, FunctionType);
90  Evaluate(const InputType & point) const override;
91 
92 protected:
94  ~BinaryThresholdSpatialFunction() override = default;
95  void
96  PrintSelf(std::ostream & os, Indent indent) const override;
97 
100 
101  typename FunctionType::Pointer m_Function{};
102 };
103 } // end namespace itk
104 
105 #ifndef ITK_MANUAL_INSTANTIATION
106 # include "itkBinaryThresholdSpatialFunction.hxx"
107 #endif
108 
109 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::BinaryThresholdSpatialFunction::FunctionOutputType
typename TFunction::OutputType FunctionOutputType
Definition: itkBinaryThresholdSpatialFunction.h:71
itk::NumericTraits::NonpositiveMin
static constexpr T NonpositiveMin()
Definition: itkNumericTraits.h:98
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::FunctionBase::OutputType
TOutput OutputType
Definition: itkFunctionBase.h:62
itk::BinaryThresholdSpatialFunction::InputType
typename TFunction::InputType InputType
Definition: itkBinaryThresholdSpatialFunction.h:65
itk::point
*par Constraints *The filter requires an image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating point(float or double) or *a user-defined "real" numerical type with arithmetic operations defined *sufficient to compute derivatives. **\par Performance *This filter will automatically multithread if run with *SetUsePrincipleComponents
itk::BinaryThresholdSpatialFunction
A spatial functions that returns if the internal spatial function is within user specified thresholds...
Definition: itkBinaryThresholdSpatialFunction.h:42
itkSpatialFunction.h
itk::NumericTraits::max
static constexpr T max(const T &)
Definition: itkNumericTraits.h:168
itk::SpatialFunction
N-dimensional spatial function class.
Definition: itkSpatialFunction.h:45
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::BinaryThresholdSpatialFunction::FunctionType
TFunction FunctionType
Definition: itkBinaryThresholdSpatialFunction.h:68
itkImageBase.h