ITK  6.0.0
Insight Toolkit
itkBinaryReconstructionLabelMapFilter.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 itkBinaryReconstructionLabelMapFilter_h
19 #define itkBinaryReconstructionLabelMapFilter_h
20 
23 
24 namespace itk
25 {
44 template <typename TImage,
45  typename TMarkerImage,
46  typename TAttributeAccessor =
47  typename Functor::AttributeLabelObjectAccessor<typename TImage::LabelObjectType>>
48 class ITK_TEMPLATE_EXPORT BinaryReconstructionLabelMapFilter : public InPlaceLabelMapFilter<TImage>
49 {
50 public:
51  ITK_DISALLOW_COPY_AND_MOVE(BinaryReconstructionLabelMapFilter);
52 
58 
60  using ImageType = TImage;
61  using ImagePointer = typename ImageType::Pointer;
63  using PixelType = typename ImageType::PixelType;
64  using IndexType = typename ImageType::IndexType;
65  using LabelObjectType = typename ImageType::LabelObjectType;
66 
67  using MarkerImageType = TMarkerImage;
70  using MarkerImagePixelType = typename MarkerImageType::PixelType;
71 
72  using AttributeAccessorType = TAttributeAccessor;
73 
75  static constexpr unsigned int ImageDimension = TImage::ImageDimension;
76 
78  itkNewMacro(Self);
79 
81  itkOverrideGetNameOfClassMacro(BinaryReconstructionLabelMapFilter);
82 
83 #ifdef ITK_USE_CONCEPT_CHECKING
84  // Begin concept checking
85  /* itkConceptMacro(InputEqualityComparableCheck,
86  (Concept::EqualityComparable<PixelType>));
87  itkConceptMacro(IntConvertibleToInputCheck,
88  (Concept::Convertible<int, PixelType>));
89  itkConceptMacro(InputOStreamWritableCheck,
90  (Concept::OStreamWritable<PixelType>));*/
91  // End concept checking
92 #endif
93 
95  itkSetInputMacro(MarkerImage, MarkerImageType);
96  itkGetInputMacro(MarkerImage, MarkerImageType);
100  void
101  SetInput1(TImage * input)
102  {
103  this->SetInput(input);
104  }
105 
107  void
108  SetInput2(TMarkerImage * input)
109  {
110  this->SetMarkerImage(input);
111  }
112 
117  itkSetMacro(ForegroundValue, MarkerImagePixelType);
118  itkGetConstMacro(ForegroundValue, MarkerImagePixelType);
121 protected:
123  ~BinaryReconstructionLabelMapFilter() override = default;
124 
125  void
126  ThreadedProcessLabelObject(LabelObjectType * labelObject) override;
127 
128  void
129  PrintSelf(std::ostream & os, Indent indent) const override;
130 
131 private:
132  MarkerImagePixelType m_ForegroundValue{};
133 
134 }; // end of class
135 
136 } // end namespace itk
137 
138 #ifndef ITK_MANUAL_INSTANTIATION
139 # include "itkBinaryReconstructionLabelMapFilter.hxx"
140 #endif
141 
142 #endif
itk::BinaryReconstructionLabelMapFilter::PixelType
typename ImageType::PixelType PixelType
Definition: itkBinaryReconstructionLabelMapFilter.h:63
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::BinaryReconstructionLabelMapFilter::AttributeAccessorType
TAttributeAccessor AttributeAccessorType
Definition: itkBinaryReconstructionLabelMapFilter.h:72
itk::InPlaceLabelMapFilter
Base class for filters that takes an image as input and overwrites that image as the output.
Definition: itkInPlaceLabelMapFilter.h:83
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::BinaryReconstructionLabelMapFilter::MarkerImageType
TMarkerImage MarkerImageType
Definition: itkBinaryReconstructionLabelMapFilter.h:67
itk::BinaryReconstructionLabelMapFilter::LabelObjectType
typename ImageType::LabelObjectType LabelObjectType
Definition: itkBinaryReconstructionLabelMapFilter.h:65
itk::BinaryReconstructionLabelMapFilter::ImagePointer
typename ImageType::Pointer ImagePointer
Definition: itkBinaryReconstructionLabelMapFilter.h:61
itk::BinaryReconstructionLabelMapFilter
Mark the objects at least partially at the same position as the objects in a binary image.
Definition: itkBinaryReconstructionLabelMapFilter.h:48
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::BinaryReconstructionLabelMapFilter::MarkerImageConstPointer
typename MarkerImageType::ConstPointer MarkerImageConstPointer
Definition: itkBinaryReconstructionLabelMapFilter.h:69
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itkAttributeLabelObject.h
itk::BinaryReconstructionLabelMapFilter::IndexType
typename ImageType::IndexType IndexType
Definition: itkBinaryReconstructionLabelMapFilter.h:64
itk::BinaryReconstructionLabelMapFilter::SetInput2
void SetInput2(TMarkerImage *input)
Definition: itkBinaryReconstructionLabelMapFilter.h:108
itk::BinaryReconstructionLabelMapFilter::MarkerImagePointer
typename MarkerImageType::Pointer MarkerImagePointer
Definition: itkBinaryReconstructionLabelMapFilter.h:68
itk::BinaryReconstructionLabelMapFilter::ImageType
TImage ImageType
Definition: itkBinaryReconstructionLabelMapFilter.h:60
itkInPlaceLabelMapFilter.h
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::BinaryReconstructionLabelMapFilter::MarkerImagePixelType
typename MarkerImageType::PixelType MarkerImagePixelType
Definition: itkBinaryReconstructionLabelMapFilter.h:70
itk::BinaryReconstructionLabelMapFilter::SetInput1
void SetInput1(TImage *input)
Definition: itkBinaryReconstructionLabelMapFilter.h:101
itk::BinaryReconstructionLabelMapFilter::ImageConstPointer
typename ImageType::ConstPointer ImageConstPointer
Definition: itkBinaryReconstructionLabelMapFilter.h:62