ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkFastMarchingImageFilter.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 itkFastMarchingImageFilter_h
19 #define itkFastMarchingImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
23 #include "itkLevelSet.h"
24 #include "itkMath.h"
25 
26 #include <functional>
27 #include <queue>
28 #include "itkMath.h"
29 
30 namespace itk
31 {
106 template<
107  typename TLevelSet,
108  typename TSpeedImage = Image< float, TLevelSet ::ImageDimension > >
109 class ITK_TEMPLATE_EXPORT FastMarchingImageFilter:
110  public ImageToImageFilter< TSpeedImage, TLevelSet >
111 {
112 public:
118 
120  itkNewMacro(Self);
121 
124 
135  typedef typename LevelSetImageType::RegionType OutputRegionType;
136  typedef typename LevelSetImageType::SpacingType OutputSpacingType;
139 
140  class AxisNodeType:public NodeType
141  {
142 public:
143  AxisNodeType() : m_Axis(0) {}
144  int GetAxis() const { return m_Axis; }
145  void SetAxis(int axis) { m_Axis = axis; }
146  const AxisNodeType & operator=(const NodeType & node)
147  { this->NodeType::operator=(node); return *this; }
148 
149 private:
150  int m_Axis;
151  };
152 
154  typedef TSpeedImage SpeedImageType;
155 
157  typedef typename SpeedImageType::Pointer SpeedImagePointer;
158  typedef typename SpeedImageType::ConstPointer SpeedImageConstPointer;
159 
161  itkStaticConstMacro(SetDimension, unsigned int,
162  LevelSetType::SetDimension);
163  itkStaticConstMacro(SpeedImageDimension, unsigned int,
164  SpeedImageType::ImageDimension);
166 
169 
174  enum LabelType { FarPoint = 0, AlivePoint,
175  TrialPoint, InitialTrialPoint, OutsidePoint };
176 
179 
182 
183  template< typename TPixel >
185  {
186  typedef Image< TPixel, SetDimension > InternalImageType;
188  InternalRegionIterator;
189  InternalRegionIterator b_it( iImage, iImage->GetLargestPossibleRegion() );
190  b_it.GoToBegin();
191 
192  TPixel zero_value = NumericTraits< TPixel >::ZeroValue();
193  typename NodeContainer::ElementIdentifier NumberOfPoints = 0;
194 
195  NodeType node;
196  node.SetValue( 0. );
197 
198  while( !b_it.IsAtEnd() )
199  {
200  if( Math::ExactlyEquals(b_it.Get(), zero_value) )
201  {
202  if( NumberOfPoints == 0 )
203  {
204  m_OutsidePoints = NodeContainer::New();
205  }
206  node.SetIndex( b_it.GetIndex() );
207  m_OutsidePoints->InsertElement( NumberOfPoints++, node );
208 
209  }
210  ++b_it;
211  }
212  this->Modified();
213  }
214 
217  {
218  m_OutsidePoints = points;
219  this->Modified();
220  }
222 
226  {
227  m_AlivePoints = points;
228  this->Modified();
229  }
231 
234  {
235  return m_AlivePoints;
236  }
237 
241  {
242  m_TrialPoints = points;
243  this->Modified();
244  }
246 
249  {
250  return m_TrialPoints;
251  }
252 
255  {
256  return m_LabelImage;
257  }
258 
262  void SetSpeedConstant(double value)
263  {
264  m_SpeedConstant = value;
265  m_InverseSpeed = -1.0 * itk::Math::sqr(1.0 / m_SpeedConstant);
266  this->Modified();
267  }
269 
271  itkGetConstReferenceMacro(SpeedConstant, double);
272 
277  itkSetMacro(NormalizationFactor, double);
278  itkGetConstMacro(NormalizationFactor, double);
280 
284  itkSetMacro(StoppingValue, double);
285 
287  itkGetConstReferenceMacro(StoppingValue, double);
288 
293  itkSetMacro(CollectPoints, bool);
294 
296  itkGetConstReferenceMacro(CollectPoints, bool);
297  itkBooleanMacro(CollectPoints);
299 
305  {
306  return m_ProcessedPoints;
307  }
308 
315  virtual void SetOutputSize(const OutputSizeType & size)
316  { m_OutputRegion = size; }
318  { return m_OutputRegion.GetSize(); }
319  itkSetMacro(OutputRegion, OutputRegionType);
320  itkGetConstReferenceMacro(OutputRegion, OutputRegionType);
321  itkSetMacro(OutputSpacing, OutputSpacingType);
322  itkGetConstReferenceMacro(OutputSpacing, OutputSpacingType);
323  itkSetMacro(OutputDirection, OutputDirectionType);
324  itkGetConstReferenceMacro(OutputDirection, OutputDirectionType);
325  itkSetMacro(OutputOrigin, OutputPointType);
326  itkGetConstReferenceMacro(OutputOrigin, OutputPointType);
327  itkSetMacro(OverrideOutputInformation, bool);
328  itkGetConstReferenceMacro(OverrideOutputInformation, bool);
329  itkBooleanMacro(OverrideOutputInformation);
331 
332 #ifdef ITK_USE_CONCEPT_CHECKING
333  // Begin concept checking
334  itkConceptMacro( SameDimensionCheck,
336  itkConceptMacro( SpeedConvertibleToDoubleCheck,
338  itkConceptMacro( DoubleConvertibleToLevelSetCheck,
340  itkConceptMacro( LevelSetOStreamWritableCheck,
342  // End concept checking
343 #endif
344 
345 protected:
347  ~FastMarchingImageFilter() ITK_OVERRIDE {}
348  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
349 
350  virtual void Initialize(LevelSetImageType *);
351 
352  virtual void UpdateNeighbors(const IndexType & index,
353  const SpeedImageType *, LevelSetImageType *);
354 
355  virtual double UpdateValue(const IndexType & index,
356  const SpeedImageType *, LevelSetImageType *);
357 
358  const AxisNodeType & GetNodeUsedInCalculation(unsigned int idx) const
359  { return m_NodesUsed[idx]; }
360 
361  void GenerateData() ITK_OVERRIDE;
362 
364  virtual void GenerateOutputInformation() ITK_OVERRIDE;
365 
366  virtual void EnlargeOutputRequestedRegion(DataObject *output) ITK_OVERRIDE;
367 
372  itkGetConstReferenceMacro(LargeValue, PixelType);
373 
374  OutputRegionType m_BufferedRegion;
376  LevelSetIndexType m_StartIndex;
377  LevelSetIndexType m_LastIndex;
378 
379  itkGetConstReferenceMacro(StartIndex, LevelSetIndexType);
380  itkGetConstReferenceMacro(LastIndex, LevelSetIndexType);
381 
382 private:
383  ITK_DISALLOW_COPY_AND_ASSIGN(FastMarchingImageFilter);
384 
385  NodeContainerPointer m_AlivePoints;
386  NodeContainerPointer m_TrialPoints;
387  NodeContainerPointer m_OutsidePoints;
388 
389  LabelImagePointer m_LabelImage;
390 
391  double m_SpeedConstant;
392  double m_InverseSpeed;
393  double m_StoppingValue;
394 
395  bool m_CollectPoints;
396  NodeContainerPointer m_ProcessedPoints;
397 
398  OutputRegionType m_OutputRegion;
399  OutputPointType m_OutputOrigin;
400  OutputSpacingType m_OutputSpacing;
401  OutputDirectionType m_OutputDirection;
402  bool m_OverrideOutputInformation;
403 
404  typename LevelSetImageType::PixelType m_LargeValue;
405  AxisNodeType m_NodesUsed[SetDimension];
406 
410  typedef std::vector< AxisNodeType > HeapContainer;
411  typedef std::greater< AxisNodeType > NodeComparer;
412  typedef std::priority_queue< AxisNodeType, HeapContainer, NodeComparer >
414 
415  HeapType m_TrialHeap;
416 
417  double m_NormalizationFactor;
418 };
419 } // namespace itk
420 
421 #ifndef ITK_MANUAL_INSTANTIATION
422 #include "itkFastMarchingImageFilter.hxx"
423 #endif
424 
425 #endif
void SetBinaryMask(Image< TPixel, SetDimension > *iImage)
SpeedImageType::ConstPointer SpeedImageConstPointer
Light weight base class for most itk classes.
TLevelSet::PixelType PixelType
Definition: itkLevelSet.h:55
std::greater< AxisNodeType > NodeComparer
TElementIdentifier ElementIdentifier
std::vector< AxisNodeType > HeapContainer
virtual const RegionType & GetLargestPossibleRegion() const
Definition: itkImageBase.h:259
void SetValue(const PixelType &input)
LevelSetImageType::RegionType OutputRegionType
LabelImageType::Pointer LabelImagePointer
virtual void SetOutputSize(const OutputSizeType &size)
LevelSetImageType::IndexType LevelSetIndexType
void SetOutsidePoints(NodeContainer *points)
std::priority_queue< AxisNodeType, HeapContainer, NodeComparer > HeapType
bool ExactlyEquals(const TInput1 &x1, const TInput2 &x2)
Return the result of an exact comparison between two scalar values of potetially different types...
Definition: itkMath.h:710
void SetTrialPoints(NodeContainer *points)
A multi-dimensional iterator templated over image type that walks an image region and is specialized ...
virtual OutputSizeType GetOutputSize() const
SmartPointer< const Self > ConstPointer
void SetIndex(const IndexType &input)
LevelSetImageType::DirectionType OutputDirectionType
LevelSetImageType::SpacingType OutputSpacingType
Image< unsigned char, itkGetStaticConstMacro(SetDimension) > LabelImageType
bool sqr(const bool x)
Definition: itkMath.h:824
Represent a node in a level set.
LevelSetTypeDefault< TLevelSet > LevelSetType
const AxisNodeType & GetNodeUsedInCalculation(unsigned int idx) const
LevelSetType::NodeContainer NodeContainer
const AxisNodeType & operator=(const NodeType &node)
LevelSetType::LevelSetPointer LevelSetPointer
LabelImagePointer GetLabelImage() const
Base class for filters that take an image as input and produce an image as output.
Define a front-end to the STL &quot;vector&quot; container that conforms to the IndexedContainerInterface.
LevelSetType::LevelSetImageType LevelSetImageType
SpeedImageType::Pointer SpeedImagePointer
Control indentation during Print() invocation.
Definition: itkIndent.h:49
LevelSetImageType::PointType OutputPointType
Solve an Eikonal equation using Fast Marching.
LevelSetImageType::SizeType OutputSizeType
#define itkConceptMacro(name, concept)
Level set type information.
Definition: itkLevelSet.h:40
Index< itkGetStaticConstMacro(SetDimension) > IndexType
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75
void SetAlivePoints(NodeContainer *points)
LevelSetType::NodeContainerPointer NodeContainerPointer
NodeContainerPointer GetProcessedPoints() const
TLevelSet::Pointer LevelSetPointer
Definition: itkLevelSet.h:51