ITK  5.4.0
Insight Toolkit
itkBinaryMaskToNarrowBandPointSetFilter.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 itkBinaryMaskToNarrowBandPointSetFilter_h
19#define itkBinaryMaskToNarrowBandPointSetFilter_h
20
25
26namespace itk
27{
53template <typename TInputImage, typename TOutputMesh>
54class ITK_TEMPLATE_EXPORT BinaryMaskToNarrowBandPointSetFilter : public ImageToMeshFilter<TInputImage, TOutputMesh>
55{
56public:
57 ITK_DISALLOW_COPY_AND_MOVE(BinaryMaskToNarrowBandPointSetFilter);
58
61
63
66
68 itkNewMacro(Self);
69
71 itkOverrideGetNameOfClassMacro(BinaryMaskToNarrowBandPointSetFilter);
72
74 using InputImageType = TInputImage;
77 using InputImagePixelType = typename InputImageType::PixelType;
78
80
82 using OutputMeshType = TOutputMesh;
86 using PointsContainer = typename OutputMeshType::PointsContainer;
87 using PointIdentifier = typename OutputMeshType::PointIdentifier;
89 using PointsContainerIterator = typename PointsContainer::Iterator;
90 using PointDataContainer = typename OutputMeshType::PointDataContainer;
92 using PointDataContainerIterator = typename PointDataContainer::Iterator;
93
95 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
96
99
108 using NodeType = typename NodeContainer::Element;
109
114
116
118 static constexpr unsigned int PointDimension = TOutputMesh::PointDimension;
119
121 using Superclass::SetInput;
122 void
123 SetInput(const InputImageType * inputImage);
124
129 itkSetMacro(BandWidth, float);
130 itkGetConstMacro(BandWidth, float);
133protected:
136 void
137 PrintSelf(std::ostream & os, Indent indent) const override;
138
140 void
141 GenerateData() override;
142
144 void
146
147private:
148 DistanceFilterPointer m_DistanceFilter{};
149 RescaleFilterPointer m_RescaleFilter{};
150
151 float m_BandWidth{};
152};
153} // end namespace itk
154
155#ifndef ITK_MANUAL_INSTANTIATION
156# include "itkBinaryMaskToNarrowBandPointSetFilter.hxx"
157#endif
158
159#endif
Generate a PointSet containing the narrow band around the edges of a input binary image.
typename DistanceFilterType::NodeContainer NodeContainer
typename PointDataContainer::Iterator PointDataContainerIterator
void SetInput(const InputImageType *inputImage)
typename OutputMeshType::PointIdentifier PointIdentifier
typename OutputMeshType::PointsContainer PointsContainer
typename DistanceFilterType::NodeContainerPointer NodeContainerPointer
typename OutputMeshType::PointDataContainer PointDataContainer
~BinaryMaskToNarrowBandPointSetFilter() override=default
void PrintSelf(std::ostream &os, Indent indent) const override
A multi-dimensional iterator templated over image type that walks an image region and is specialized ...
ImageToMeshFilter is the base class for all process objects that output Mesh data and require image d...
typename InputImageType::ConstPointer InputImageConstPointer
typename InputImageType::RegionType InputImageRegionType
typename InputImageType::PixelType InputImagePixelType
typename OutputMeshType::Pointer OutputMeshPointer
Templated n-dimensional image class.
Definition: itkImage.h:89
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
Reinitialize the level set to the signed distance function.
typename LevelSetType::NodeContainer NodeContainer
typename LevelSetType::NodeContainerPointer NodeContainerPointer
Applies a linear transformation to the intensity levels of the input Image.
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....