ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkBinaryMaskToNarrowBandPointSetFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software 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 #ifndef itkBinaryMaskToNarrowBandPointSetFilter_h
19 #define itkBinaryMaskToNarrowBandPointSetFilter_h
20 
21 #include "itkImageToMeshFilter.h"
25 
26 namespace itk
27 {
53 template< typename TInputImage, typename TOutputMesh >
54 class ITK_TEMPLATE_EXPORT BinaryMaskToNarrowBandPointSetFilter:
55  public ImageToMeshFilter< TInputImage, TOutputMesh >
56 {
57 public:
60 
62 
65 
67  itkNewMacro(Self);
68 
71 
73  typedef TInputImage InputImageType;
74  typedef typename InputImageType::ConstPointer InputImageConstPointer;
75  typedef typename InputImageType::RegionType InputImageRegionType;
76  typedef typename InputImageType::PixelType InputImagePixelType;
77 
79 
81  typedef TOutputMesh OutputMeshType;
83  typedef typename OutputMeshType::Pointer OutputMeshPointer;
84  typedef typename OutputMeshType::ConstPointer OutputMeshConstPointer;
85  typedef typename OutputMeshType::PointsContainer PointsContainer;
86  typedef typename OutputMeshType::PointIdentifier PointIdentifier;
87  typedef typename PointsContainer::Pointer PointsContainerPointer;
88  typedef typename PointsContainer::Iterator PointsContainerIterator;
89  typedef typename OutputMeshType::PointDataContainer PointDataContainer;
90  typedef typename PointDataContainer::Pointer PointDataContainerPointer;
91  typedef typename PointDataContainer::Iterator PointDataContainerIterator;
92 
94  itkStaticConstMacro(ImageDimension, unsigned int,
95  TInputImage::ImageDimension);
96 
98  typedef itk::Image< float,
99  itkGetStaticConstMacro(ImageDimension) > RealImageType;
100 
110 
115 
117 
119  itkStaticConstMacro(PointDimension, unsigned int,
120  TOutputMesh::PointDimension);
121 
123  void GenerateData(void) ITK_OVERRIDE;
124 
126  void GenerateOutputInformation(void) ITK_OVERRIDE;
127 
129  using Superclass::SetInput;
130  void SetInput(const InputImageType *inputImage);
131 
136  itkSetMacro(BandWidth, float);
137  itkGetConstMacro(BandWidth, float);
139 
140 protected:
142  ~BinaryMaskToNarrowBandPointSetFilter() ITK_OVERRIDE;
143  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
144 
145 private:
146  ITK_DISALLOW_COPY_AND_ASSIGN(BinaryMaskToNarrowBandPointSetFilter);
147 
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
Light weight base class for most itk classes.
DistanceFilterType::NodeContainerPointer NodeContainerPointer
ImageToMeshFilter< TInputImage, TOutputMesh > Superclass
ImageRegionConstIteratorWithIndex< InputImageType > InputImageIterator
ReinitializeLevelSetImageFilter< RealImageType > DistanceFilterType
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...
Generate a PointSet containing the narrow band around the edges of a input binary image...
RescaleIntensityImageFilter< InputImageType, RealImageType > RescaleFilterType
itk::Image< float, itkGetStaticConstMacro(ImageDimension) > RealImageType
Reinitialize the level set to the signed distance function.
Applies a linear transformation to the intensity levels of the input Image.
Define a front-end to the STL &quot;vector&quot; container that conforms to the IndexedContainerInterface.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Templated n-dimensional image class.
Definition: itkImage.h:75