ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkFastMarchingUpwindGradientImageFilter.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 itkFastMarchingUpwindGradientImageFilter_h
19 #define itkFastMarchingUpwindGradientImageFilter_h
20 
22 #include "itkImage.h"
23 
24 namespace itk
25 {
59 template<
60  typename TLevelSet,
61  typename TSpeedImage = Image< float, TLevelSet ::ImageDimension > >
62 class ITK_TEMPLATE_EXPORT FastMarchingUpwindGradientImageFilter:
63  public FastMarchingImageFilter< TLevelSet, TSpeedImage >
64 {
65 public:
71 
73  itkNewMacro(Self);
74 
77 
79  typedef typename Superclass::LevelSetType LevelSetType;
80  typedef typename Superclass::SpeedImageType SpeedImageType;
81  typedef typename Superclass::LevelSetImageType LevelSetImageType;
82  typedef typename Superclass::LevelSetPointer LevelSetPointer;
83  typedef typename Superclass::SpeedImageConstPointer SpeedImageConstPointer;
84  typedef typename Superclass::LabelImageType LabelImageType;
85  typedef typename Superclass::PixelType PixelType;
86  typedef typename Superclass::AxisNodeType AxisNodeType;
87  typedef typename Superclass::NodeType NodeType;
88  typedef typename Superclass::NodeContainer NodeContainer;
90 
91  typedef typename Superclass::IndexType IndexType;
92  typedef typename Superclass::OutputSpacingType OutputSpacingType;
93  typedef typename Superclass::LevelSetIndexType LevelSetIndexType;
94 
95  typedef typename Superclass::OutputPointType PointType;
96 
98  itkStaticConstMacro(SetDimension, unsigned int, Superclass::SetDimension);
99 
104  {
105  m_TargetPoints = points;
106  this->Modified();
107  }
109 
112  { return m_TargetPoints; }
113 
116  { return m_ReachedTargetPoints; }
117 
119  typedef CovariantVector< PixelType,
120  itkGetStaticConstMacro(SetDimension) > GradientPixelType;
121 
123  typedef Image< GradientPixelType,
124  itkGetStaticConstMacro(SetDimension) > GradientImageType;
125 
128 
131  { return m_GradientImage; }
132 
135  itkSetMacro(GenerateGradientImage, bool);
136 
138  itkGetConstReferenceMacro(GenerateGradientImage, bool);
139  itkBooleanMacro(GenerateGradientImage);
141 
145  itkSetMacro(TargetOffset, double);
146 
148  itkGetConstReferenceMacro(TargetOffset, double);
149 
153  itkSetMacro(TargetReachedMode, int);
154  itkGetConstReferenceMacro(TargetReachedMode, int);
156  { this->SetTargetReachedMode(NoTargets); }
158  { this->SetTargetReachedMode(OneTarget); }
160  {
161  this->SetTargetReachedMode(SomeTargets);
162  m_NumberOfTargets = numberOfTargets;
163  }
165 
167  { this->SetTargetReachedMode(AllTargets); }
168 
170  itkGetConstReferenceMacro(NumberOfTargets, SizeValueType);
171 
176  itkGetConstReferenceMacro(TargetValue, double);
177 
178  enum {
182  AllTargets
183  };
184 
185 #ifdef ITK_USE_CONCEPT_CHECKING
186  // Begin concept checking
187  itkConceptMacro( LevelSetDoubleDivisionOperatorsCheck,
189  itkConceptMacro( LevelSetDoubleDivisionAndAssignOperatorsCheck,
191  // End concept checking
192 #endif
193 
194 protected:
197  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
198 
199  virtual void Initialize(LevelSetImageType *) ITK_OVERRIDE;
200 
201  void GenerateData() ITK_OVERRIDE;
202 
203  virtual void UpdateNeighbors(const IndexType & index,
204  const SpeedImageType *, LevelSetImageType *) ITK_OVERRIDE;
205 
206  virtual void ComputeGradient(const IndexType & index,
207  const LevelSetImageType *output,
208  const LabelImageType *labelImage,
209  GradientImageType *gradientImage);
210 
211 private:
212  ITK_DISALLOW_COPY_AND_ASSIGN(FastMarchingUpwindGradientImageFilter);
213 
214  NodeContainerPointer m_TargetPoints;
215  NodeContainerPointer m_ReachedTargetPoints;
216 
217  GradientImagePointer m_GradientImage;
218 
219  bool m_GenerateGradientImage;
220 
221  double m_TargetOffset;
222 
223  int m_TargetReachedMode;
224 
225  double m_TargetValue;
226 
227  SizeValueType m_NumberOfTargets;
228 };
229 } // namespace itk
230 
231 #ifndef ITK_MANUAL_INSTANTIATION
232 #include "itkFastMarchingUpwindGradientImageFilter.hxx"
233 #endif
234 
235 #endif
FastMarchingImageFilter< TLevelSet, TSpeedImage > Superclass
Light weight base class for most itk classes.
unsigned long SizeValueType
Definition: itkIntTypes.h:143
Image< GradientPixelType, itkGetStaticConstMacro(SetDimension) > GradientImageType
Generates the upwind gradient field of fast marching arrival times.
CovariantVector< PixelType, itkGetStaticConstMacro(SetDimension) > GradientPixelType
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
Solve an Eikonal equation using Fast Marching.
#define itkConceptMacro(name, concept)
A templated class holding a n-Dimensional covariant vector.
Templated n-dimensional image class.
Definition: itkImage.h:75