ITK  5.0.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:
66  ITK_DISALLOW_COPY_AND_ASSIGN(FastMarchingUpwindGradientImageFilter);
67 
73 
75  itkNewMacro(Self);
76 
79 
81  using LevelSetType = typename Superclass::LevelSetType;
82  using SpeedImageType = typename Superclass::SpeedImageType;
83  using LevelSetImageType = typename Superclass::LevelSetImageType;
84  using LevelSetPointer = typename Superclass::LevelSetPointer;
85  using SpeedImageConstPointer = typename Superclass::SpeedImageConstPointer;
86  using LabelImageType = typename Superclass::LabelImageType;
87  using PixelType = typename Superclass::PixelType;
88  using AxisNodeType = typename Superclass::AxisNodeType;
89  using NodeType = typename Superclass::NodeType;
90  using NodeContainer = typename Superclass::NodeContainer;
92 
93  using IndexType = typename Superclass::IndexType;
94  using OutputSpacingType = typename Superclass::OutputSpacingType;
95  using LevelSetIndexType = typename Superclass::LevelSetIndexType;
96 
97  using PointType = typename Superclass::OutputPointType;
98 
100  static constexpr unsigned int SetDimension = Superclass::SetDimension;
101 
106  {
107  m_TargetPoints = points;
108  this->Modified();
109  }
111 
114  { return m_TargetPoints; }
115 
118  { return m_ReachedTargetPoints; }
119 
121  using GradientPixelType = CovariantVector< PixelType,
122  Self::SetDimension >;
123 
126  Self::SetDimension >;
127 
130 
133  { return m_GradientImage; }
134 
137  itkSetMacro(GenerateGradientImage, bool);
138 
140  itkGetConstReferenceMacro(GenerateGradientImage, bool);
141  itkBooleanMacro(GenerateGradientImage);
143 
147  itkSetMacro(TargetOffset, double);
148 
150  itkGetConstReferenceMacro(TargetOffset, double);
151 
155  itkSetMacro(TargetReachedMode, int);
156  itkGetConstReferenceMacro(TargetReachedMode, int);
158  { this->SetTargetReachedMode(NoTargets); }
160  { this->SetTargetReachedMode(OneTarget); }
162  {
163  this->SetTargetReachedMode(SomeTargets);
164  m_NumberOfTargets = numberOfTargets;
165  }
167 
169  { this->SetTargetReachedMode(AllTargets); }
170 
172  itkGetConstReferenceMacro(NumberOfTargets, SizeValueType);
173 
178  itkGetConstReferenceMacro(TargetValue, double);
179 
180  enum {
184  AllTargets
185  };
186 
187 #ifdef ITK_USE_CONCEPT_CHECKING
188  // Begin concept checking
189  itkConceptMacro( LevelSetDoubleDivisionOperatorsCheck,
191  itkConceptMacro( LevelSetDoubleDivisionAndAssignOperatorsCheck,
193  // End concept checking
194 #endif
195 
196 protected:
198  ~FastMarchingUpwindGradientImageFilter() override = default;
199  void PrintSelf(std::ostream & os, Indent indent) const override;
200 
201  void Initialize(LevelSetImageType *) override;
202 
203  void GenerateData() override;
204 
205  void UpdateNeighbors(const IndexType & index,
206  const SpeedImageType *, LevelSetImageType *) override;
207 
208  virtual void ComputeGradient(const IndexType & index,
209  const LevelSetImageType *output,
210  const LabelImageType *labelImage,
211  GradientImageType *gradientImage);
212 
213 private:
216 
218 
220 
222 
224 
226 
228 };
229 } // namespace itk
230 
231 #ifndef ITK_MANUAL_INSTANTIATION
232 #include "itkFastMarchingUpwindGradientImageFilter.hxx"
233 #endif
234 
235 #endif
typename LevelSetType::NodeType NodeType
Light weight base class for most itk classes.
typename LevelSetType::LevelSetPointer LevelSetPointer
typename LevelSetType::LevelSetImageType LevelSetImageType
Represent a n-dimensional index in a n-dimensional image.
Definition: itkIndex.h:66
unsigned long SizeValueType
Definition: itkIntTypes.h:83
typename LevelSetType::NodeContainerPointer NodeContainerPointer
Generates the upwind gradient field of fast marching arrival times.
typename LevelSetType::NodeContainer NodeContainer
typename LevelSetType::PixelType PixelType
typename SpeedImageType::ConstPointer SpeedImageConstPointer
Control indentation during Print() invocation.
Definition: itkIndent.h:49
typename LevelSetImageType::SpacingType OutputSpacingType
Solve an Eikonal equation using Fast Marching.
typename LevelSetImageType::IndexType LevelSetIndexType
#define itkConceptMacro(name, concept)
Level set type information.
Definition: itkLevelSet.h:40
A templated class holding a n-Dimensional covariant vector.
Templated n-dimensional image class.
Definition: itkImage.h:75