ITK  4.11.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 {
58 template<
59  typename TLevelSet,
60  typename TSpeedImage = Image< float, TLevelSet ::ImageDimension > >
61 class ITK_TEMPLATE_EXPORT FastMarchingUpwindGradientImageFilter:
62  public FastMarchingImageFilter< TLevelSet, TSpeedImage >
63 {
64 public:
70 
72  itkNewMacro(Self);
73 
76 
78  typedef typename Superclass::LevelSetType LevelSetType;
79  typedef typename Superclass::SpeedImageType SpeedImageType;
80  typedef typename Superclass::LevelSetImageType LevelSetImageType;
81  typedef typename Superclass::LevelSetPointer LevelSetPointer;
82  typedef typename Superclass::SpeedImageConstPointer SpeedImageConstPointer;
83  typedef typename Superclass::LabelImageType LabelImageType;
84  typedef typename Superclass::PixelType PixelType;
85  typedef typename Superclass::AxisNodeType AxisNodeType;
86  typedef typename Superclass::NodeType NodeType;
87  typedef typename Superclass::NodeContainer NodeContainer;
89 
90  typedef typename Superclass::IndexType IndexType;
91  typedef typename Superclass::OutputSpacingType OutputSpacingType;
92  typedef typename Superclass::LevelSetIndexType LevelSetIndexType;
93 
94  typedef typename Superclass::OutputPointType PointType;
95 
97  itkStaticConstMacro(SetDimension, unsigned int, Superclass::SetDimension);
98 
103  {
104  m_TargetPoints = points;
105  this->Modified();
106  }
108 
111  { return m_TargetPoints; }
112 
115  { return m_ReachedTargetPoints; }
116 
118  typedef CovariantVector< PixelType,
119  itkGetStaticConstMacro(SetDimension) > GradientPixelType;
120 
122  typedef Image< GradientPixelType,
123  itkGetStaticConstMacro(SetDimension) > GradientImageType;
124 
127 
130  { return m_GradientImage; }
131 
134  itkSetMacro(GenerateGradientImage, bool);
135 
137  itkGetConstReferenceMacro(GenerateGradientImage, bool);
138  itkBooleanMacro(GenerateGradientImage);
140 
144  itkSetMacro(TargetOffset, double);
145 
147  itkGetConstReferenceMacro(TargetOffset, double);
148 
152  itkSetMacro(TargetReachedMode, int);
153  itkGetConstReferenceMacro(TargetReachedMode, int);
155  { this->SetTargetReachedMode(NoTargets); }
157  { this->SetTargetReachedMode(OneTarget); }
159  {
160  this->SetTargetReachedMode(SomeTargets);
161  m_NumberOfTargets = numberOfTargets;
162  }
164 
166  { this->SetTargetReachedMode(AllTargets); }
167 
169  itkGetConstReferenceMacro(NumberOfTargets, SizeValueType);
170 
175  itkGetConstReferenceMacro(TargetValue, double);
176 
177  enum {
181  AllTargets
182  };
183 
184 #ifdef ITK_USE_CONCEPT_CHECKING
185  // Begin concept checking
186  itkConceptMacro( LevelSetDoubleDivisionOperatorsCheck,
188  itkConceptMacro( LevelSetDoubleDivisionAndAssignOperatorsCheck,
190  // End concept checking
191 #endif
192 
193 protected:
196  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
197 
198  virtual void Initialize(LevelSetImageType *) ITK_OVERRIDE;
199 
200  void GenerateData() ITK_OVERRIDE;
201 
202  virtual void UpdateNeighbors(const IndexType & index,
203  const SpeedImageType *, LevelSetImageType *) ITK_OVERRIDE;
204 
205  virtual void ComputeGradient(const IndexType & index,
206  const LevelSetImageType *output,
207  const LabelImageType *labelImage,
208  GradientImageType *gradientImage);
209 
210 private:
211  ITK_DISALLOW_COPY_AND_ASSIGN(FastMarchingUpwindGradientImageFilter);
212 
213  NodeContainerPointer m_TargetPoints;
214  NodeContainerPointer m_ReachedTargetPoints;
215 
216  GradientImagePointer m_GradientImage;
217 
218  bool m_GenerateGradientImage;
219 
220  double m_TargetOffset;
221 
222  int m_TargetReachedMode;
223 
224  double m_TargetValue;
225 
226  SizeValueType m_NumberOfTargets;
227 };
228 } // namespace itk
229 
230 #ifndef ITK_MANUAL_INSTANTIATION
231 #include "itkFastMarchingUpwindGradientImageFilter.hxx"
232 #endif
233 
234 #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