ITK  4.6.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 
26 namespace itk
27 {
56 template< typename TInputImage, typename TOutputImage >
58  public ImageToImageFilter< TInputImage, TOutputImage >
59 {
60 public:
66 
68  itkNewMacro(Self);
69 
72 
76 
79  itkStaticConstMacro(ImageDimension, unsigned int,
80  TInputImage::ImageDimension);
81  itkStaticConstMacro(OutputImageDimension, unsigned int,
82  TOutputImage::ImageDimension);
84 
87  typedef typename OutputImageType::PixelType PixelType;
88  typedef typename InputImageType::PixelType InputPixelType;
89  typedef typename OutputImageType::RegionType OutputImageRegionType;
90 
91  typedef typename InputImageType::SizeType InputSizeType;
92  typedef typename OutputImageType::SizeType SizeType;
93 
94  typedef typename InputImageType::IndexType InputIndexType;
95  typedef typename OutputImageType::IndexType IndexType;
96 
97  typedef typename InputImageType::SpacingType InputSpacingType;
98 
106 
109  itkSetMacro(LevelSetValue, InputPixelType);
110  itkGetConstMacro(LevelSetValue, InputPixelType);
112 
115  itkSetMacro(FarValue, PixelType);
116  itkGetConstMacro(FarValue, PixelType);
118 
121  itkSetMacro(NarrowBanding, bool);
122  itkGetConstMacro(NarrowBanding, bool);
123  itkBooleanMacro(NarrowBanding);
125 
127  void SetNarrowBand(NarrowBandType *ptr);
128 
130  { return m_NarrowBand; }
131 
132 #ifdef ITK_USE_CONCEPT_CHECKING
133  // Begin concept checking
134  itkConceptMacro( InputEqualityComparableCheck,
136  itkConceptMacro( OutputEqualityComparableCheck,
138  itkConceptMacro( SameDimensionCheck,
140  itkConceptMacro( DoubleConvertibleToOutputCheck,
142  itkConceptMacro( InputConvertibleToOutputCheck,
144  itkConceptMacro( OutputAdditiveOperatorsCheck,
146  itkConceptMacro( InputOStreamWritableCheck,
148  itkConceptMacro( OutputOStreamWritableCheck,
150  // End concept checking
151 #endif
152 
153 protected:
156  void PrintSelf(std::ostream & os, Indent indent) const;
157 
158  void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
159  ThreadIdType threadId);
160 
161  void ThreadedGenerateDataFull(const OutputImageRegionType & outputRegionForThread,
162  ThreadIdType threadId);
163 
164  void ThreadedGenerateDataBand(const OutputImageRegionType & outputRegionForThread,
165  ThreadIdType threadId);
166 
168 
169  virtual void GenerateInputRequestedRegion();
170 
172 
175 
176  void ComputeValue( const InputNeighbordIteratorType& inNeigIt,
178  unsigned int center,
179  const std::vector< OffsetValueType >& stride );
180 
181 private:
182  IsoContourDistanceImageFilter(const Self &); //purposely not implemented
183  void operator=(const Self &); //purposely not implemented
184 
187 
189 
192  std::vector< RegionType > m_NarrowBandRegion;
193 
196 };
197 } // namespace itk
198 
199 #ifndef ITK_MANUAL_INSTANTIATION
200 #include "itkIsoContourDistanceImageFilter.hxx"
201 #endif
202 
203 #endif
void SetNarrowBand(NarrowBandType *ptr)
ImageToImageFilter< TInputImage, TOutputImage > Superclass
ConstNeighborhoodIterator< InputImageType > InputNeighbordIteratorType
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...
void PrintSelf(std::ostream &os, Indent indent) const
void ThreadedGenerateDataFull(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
NodeContainerType::iterator Iterator
Definition: itkNarrowBand.h:69
void ComputeValue(const InputNeighbordIteratorType &inNeigIt, OutputNeighborhoodIteratorType &outNeigIt, unsigned int center, const std::vector< OffsetValueType > &stride)
NodeContainerType::const_iterator ConstIterator
Definition: itkNarrowBand.h:68
Narrow Band class.
Definition: itkNarrowBand.h:51
NeighborhoodIterator< OutputImageType > OutputNeighborhoodIteratorType
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
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
TOutputImage OutputImageType
virtual void EnlargeOutputRequestedRegion(DataObject *)
#define itkConceptMacro(name, concept)
virtual void GenerateInputRequestedRegion()
void ThreadedGenerateDataBand(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
Base class for all data objects in ITK.
Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.
unsigned int ThreadIdType
Definition: itkIntTypes.h:159