ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkIsoContourDistanceImageFilter.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 itkIsoContourDistanceImageFilter_h
19 #define itkIsoContourDistanceImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include "itkNarrowBand.h"
24 #include "itkBarrier.h"
25 #include "itkNumericTraits.h"
26 
27 namespace itk
28 {
57 template< typename TInputImage, typename TOutputImage >
58 class ITK_TEMPLATE_EXPORT IsoContourDistanceImageFilter:
59  public ImageToImageFilter< TInputImage, TOutputImage >
60 {
61 public:
67 
69  itkNewMacro(Self);
70 
73 
76  typedef typename Superclass::OutputImageType OutputImageType;
77 
80  itkStaticConstMacro(ImageDimension, unsigned int,
81  TInputImage::ImageDimension);
82  itkStaticConstMacro(OutputImageDimension, unsigned int,
83  TOutputImage::ImageDimension);
85 
88  typedef typename OutputImageType::PixelType PixelType;
89  typedef typename InputImageType::PixelType InputPixelType;
91 
92  typedef typename OutputImageType::RegionType OutputImageRegionType;
93 
96 
99 
100  typedef typename InputImageType::SpacingType InputSpacingType;
101 
109 
112  itkSetMacro(LevelSetValue, PixelRealType);
113  itkGetConstMacro(LevelSetValue, PixelRealType);
115 
118  itkSetMacro(FarValue, PixelType);
119  itkGetConstMacro(FarValue, PixelType);
121 
124  itkSetMacro(NarrowBanding, bool);
125  itkGetConstMacro(NarrowBanding, bool);
126  itkBooleanMacro(NarrowBanding);
128 
130  void SetNarrowBand(NarrowBandType *ptr);
131 
133  { return m_NarrowBand; }
134 
135 #ifdef ITK_USE_CONCEPT_CHECKING
136  // Begin concept checking
137  itkConceptMacro( InputEqualityComparableCheck,
139  itkConceptMacro( OutputEqualityComparableCheck,
141  itkConceptMacro( SameDimensionCheck,
143  itkConceptMacro( DoubleConvertibleToOutputCheck,
145  itkConceptMacro( InputConvertibleToOutputCheck,
147  itkConceptMacro( OutputAdditiveOperatorsCheck,
149  itkConceptMacro( InputOStreamWritableCheck,
151  itkConceptMacro( OutputOStreamWritableCheck,
153  // End concept checking
154 #endif
155 
156 protected:
159  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
160 
161  void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
162  ThreadIdType threadId) ITK_OVERRIDE;
163 
164  void ThreadedGenerateDataFull(const OutputImageRegionType & outputRegionForThread,
165  ThreadIdType threadId);
166 
167  void ThreadedGenerateDataBand(const OutputImageRegionType & outputRegionForThread,
168  ThreadIdType threadId);
169 
170  void BeforeThreadedGenerateData() ITK_OVERRIDE;
171 
172  virtual void GenerateInputRequestedRegion() ITK_OVERRIDE;
173 
174  virtual void EnlargeOutputRequestedRegion(DataObject *) ITK_OVERRIDE;
175 
178 
179  void ComputeValue( const InputNeighbordIteratorType& inNeigIt,
180  OutputNeighborhoodIteratorType& outNeigIt,
181  unsigned int center,
182  const std::vector< OffsetValueType >& stride );
183 
184 private:
185  ITK_DISALLOW_COPY_AND_ASSIGN(IsoContourDistanceImageFilter);
186 
187  PixelRealType m_LevelSetValue;
188  PixelType m_FarValue;
189 
190  InputSpacingType m_Spacing;
191 
192  bool m_NarrowBanding;
193  NarrowBandPointer m_NarrowBand;
194  std::vector< RegionType > m_NarrowBandRegion;
195 
197  typename Barrier::Pointer m_Barrier;
198 };
199 } // namespace itk
200 
201 #ifndef ITK_MANUAL_INSTANTIATION
202 #include "itkIsoContourDistanceImageFilter.hxx"
203 #endif
204 
205 #endif
signed long OffsetValueType
Definition: itkIntTypes.h:154
ImageToImageFilter< TInputImage, TOutputImage > Superclass
Compute an approximate distance from an interpolated isocontour to the close grid points...
Base class for all process objects that output image data.
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
NodeContainerType::iterator Iterator
Definition: itkNarrowBand.h:69
NumericTraits< InputPixelType >::RealType PixelRealType
NodeContainerType::const_iterator ConstIterator
Definition: itkNarrowBand.h:68
Narrow Band class.
Definition: itkNarrowBand.h:51
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Define additional traits for native types such as int or float.
TOutputImage OutputImageType
#define itkConceptMacro(name, concept)
Base class for all data objects in ITK.
Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.
Standard barrier class implementation for synchronizing the execution of threads. ...
Definition: itkBarrier.h:41