ITK  5.1.0
Insight 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 <typename TLevelSet, typename TSpeedImage = Image<float, TLevelSet ::ImageDimension>>
60 class ITK_TEMPLATE_EXPORT FastMarchingUpwindGradientImageFilter : public FastMarchingImageFilter<TLevelSet, TSpeedImage>
61 {
62 public:
63  ITK_DISALLOW_COPY_AND_ASSIGN(FastMarchingUpwindGradientImageFilter);
64 
70 
72  itkNewMacro(Self);
73 
76 
78  using LevelSetType = typename Superclass::LevelSetType;
79  using SpeedImageType = typename Superclass::SpeedImageType;
80  using LevelSetImageType = typename Superclass::LevelSetImageType;
81  using LevelSetPointer = typename Superclass::LevelSetPointer;
82  using SpeedImageConstPointer = typename Superclass::SpeedImageConstPointer;
83  using LabelImageType = typename Superclass::LabelImageType;
84  using PixelType = typename Superclass::PixelType;
85  using AxisNodeType = typename Superclass::AxisNodeType;
86  using NodeType = typename Superclass::NodeType;
87  using NodeContainer = typename Superclass::NodeContainer;
89 
90  using IndexType = typename Superclass::IndexType;
91  using OutputSpacingType = typename Superclass::OutputSpacingType;
92  using LevelSetIndexType = typename Superclass::LevelSetIndexType;
93 
94  using PointType = typename Superclass::OutputPointType;
95 
97  static constexpr unsigned int SetDimension = Superclass::SetDimension;
98 
102  void
104  {
105  m_TargetPoints = points;
106  this->Modified();
107  }
109 
111  NodeContainerPointer
113  {
114  return m_TargetPoints;
115  }
116 
118  NodeContainerPointer
120  {
121  return m_ReachedTargetPoints;
122  }
123 
126 
129 
132 
136  {
137  return m_GradientImage;
138  }
139 
142  itkSetMacro(GenerateGradientImage, bool);
143 
145  itkGetConstReferenceMacro(GenerateGradientImage, bool);
146  itkBooleanMacro(GenerateGradientImage);
148 
152  itkSetMacro(TargetOffset, double);
153 
155  itkGetConstReferenceMacro(TargetOffset, double);
156 
160  itkSetMacro(TargetReachedMode, int);
161  itkGetConstReferenceMacro(TargetReachedMode, int);
162  void
164  {
165  this->SetTargetReachedMode(NoTargets);
166  }
167  void
169  {
170  this->SetTargetReachedMode(OneTarget);
171  }
172  void
174  {
175  this->SetTargetReachedMode(SomeTargets);
176  m_NumberOfTargets = numberOfTargets;
177  }
179 
180  void
182  {
183  this->SetTargetReachedMode(AllTargets);
184  }
185 
187  itkGetConstReferenceMacro(NumberOfTargets, SizeValueType);
188 
193  itkGetConstReferenceMacro(TargetValue, double);
194 
195  enum
196  {
200  AllTargets
201  };
202 
203 #ifdef ITK_USE_CONCEPT_CHECKING
204  // Begin concept checking
205  itkConceptMacro(LevelSetDoubleDivisionOperatorsCheck,
207  itkConceptMacro(LevelSetDoubleDivisionAndAssignOperatorsCheck,
209  // End concept checking
210 #endif
211 
212 protected:
214  ~FastMarchingUpwindGradientImageFilter() override = default;
215  void
216  PrintSelf(std::ostream & os, Indent indent) const override;
217 
218  void
219  Initialize(LevelSetImageType *) override;
220 
221  void
222  GenerateData() override;
223 
224  void
225  UpdateNeighbors(const IndexType & index, const SpeedImageType *, LevelSetImageType *) override;
226 
227  virtual void
228  ComputeGradient(const IndexType & index,
229  const LevelSetImageType * output,
230  const LabelImageType * labelImage,
231  GradientImageType * gradientImage);
232 
233 private:
236 
238 
240 
242 
244 
246 
248 };
249 } // namespace itk
250 
251 #ifndef ITK_MANUAL_INSTANTIATION
252 # include "itkFastMarchingUpwindGradientImageFilter.hxx"
253 #endif
254 
255 #endif
itk::FastMarchingUpwindGradientImageFilter::SomeTargets
Definition: itkFastMarchingUpwindGradientImageFilter.h:199
itk::FastMarchingImageFilter::PixelType
typename LevelSetType::PixelType PixelType
Definition: itkFastMarchingImageFilter.h:128
itk::LevelSetTypeDefault
Level set type information.
Definition: itkLevelSet.h:40
itk::Index
Represent a n-dimensional index in a n-dimensional image.
Definition: itkIndex.h:66
itk::FastMarchingUpwindGradientImageFilter::m_TargetReachedMode
int m_TargetReachedMode
Definition: itkFastMarchingUpwindGradientImageFilter.h:243
itk::FastMarchingUpwindGradientImageFilter::SetTargetReachedModeToAllTargets
void SetTargetReachedModeToAllTargets()
Definition: itkFastMarchingUpwindGradientImageFilter.h:181
itk::FastMarchingUpwindGradientImageFilter::m_GenerateGradientImage
bool m_GenerateGradientImage
Definition: itkFastMarchingUpwindGradientImageFilter.h:239
itk::FastMarchingImageFilter::LevelSetPointer
typename LevelSetType::LevelSetPointer LevelSetPointer
Definition: itkFastMarchingImageFilter.h:127
itk::FastMarchingUpwindGradientImageFilter::m_GradientImage
GradientImagePointer m_GradientImage
Definition: itkFastMarchingUpwindGradientImageFilter.h:237
itk::FastMarchingImageFilter::OutputSpacingType
typename LevelSetImageType::SpacingType OutputSpacingType
Definition: itkFastMarchingImageFilter.h:135
itk::FastMarchingUpwindGradientImageFilter::OneTarget
Definition: itkFastMarchingUpwindGradientImageFilter.h:198
itkImage.h
itk::FastMarchingUpwindGradientImageFilter::m_ReachedTargetPoints
NodeContainerPointer m_ReachedTargetPoints
Definition: itkFastMarchingUpwindGradientImageFilter.h:235
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::FastMarchingUpwindGradientImageFilter::GetReachedTargetPoints
NodeContainerPointer GetReachedTargetPoints()
Definition: itkFastMarchingUpwindGradientImageFilter.h:119
itk::FastMarchingUpwindGradientImageFilter
Generates the upwind gradient field of fast marching arrival times.
Definition: itkFastMarchingUpwindGradientImageFilter.h:60
itk::Concept::DivisionAndAssignOperators
Definition: itkConceptChecking.h:500
itk::FastMarchingUpwindGradientImageFilter::GetTargetPoints
NodeContainerPointer GetTargetPoints()
Definition: itkFastMarchingUpwindGradientImageFilter.h:112
itk::FastMarchingUpwindGradientImageFilter::m_NumberOfTargets
SizeValueType m_NumberOfTargets
Definition: itkFastMarchingUpwindGradientImageFilter.h:247
itk::FastMarchingImageFilter::LevelSetIndexType
typename LevelSetImageType::IndexType LevelSetIndexType
Definition: itkFastMarchingImageFilter.h:406
itk::FastMarchingUpwindGradientImageFilter::SetTargetReachedModeToNoTargets
void SetTargetReachedModeToNoTargets()
Definition: itkFastMarchingUpwindGradientImageFilter.h:163
itkFastMarchingImageFilter.h
itk::FastMarchingImageFilter::NodeType
typename LevelSetType::NodeType NodeType
Definition: itkFastMarchingImageFilter.h:129
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::FastMarchingUpwindGradientImageFilter::GetGradientImage
GradientImagePointer GetGradientImage() const
Definition: itkFastMarchingUpwindGradientImageFilter.h:135
itk::FastMarchingUpwindGradientImageFilter::SetTargetPoints
void SetTargetPoints(NodeContainer *points)
Definition: itkFastMarchingUpwindGradientImageFilter.h:103
itk::FastMarchingUpwindGradientImageFilter::NoTargets
Definition: itkFastMarchingUpwindGradientImageFilter.h:197
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:59
itk::FastMarchingImageFilter
Solve an Eikonal equation using Fast Marching.
Definition: itkFastMarchingImageFilter.h:107
itk::FastMarchingImageFilter::SpeedImageType
TSpeedImage SpeedImageType
Definition: itkFastMarchingImageFilter.h:165
itk::FastMarchingUpwindGradientImageFilter::SetTargetReachedModeToSomeTargets
void SetTargetReachedModeToSomeTargets(SizeValueType numberOfTargets)
Definition: itkFastMarchingUpwindGradientImageFilter.h:173
itk::FastMarchingUpwindGradientImageFilter::m_TargetValue
double m_TargetValue
Definition: itkFastMarchingUpwindGradientImageFilter.h:245
itk::Concept::DivisionOperators
Definition: itkConceptChecking.h:471
itk::FastMarchingUpwindGradientImageFilter::PointType
typename Superclass::OutputPointType PointType
Definition: itkFastMarchingUpwindGradientImageFilter.h:94
itk::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition: itkCovariantVector.h:68
itk::FastMarchingUpwindGradientImageFilter::SetTargetReachedModeToOneTarget
void SetTargetReachedModeToOneTarget()
Definition: itkFastMarchingUpwindGradientImageFilter.h:168
itk::FastMarchingUpwindGradientImageFilter::GradientImagePointer
typename GradientImageType::Pointer GradientImagePointer
Definition: itkFastMarchingUpwindGradientImageFilter.h:131
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:64
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkArray.h:26
itk::FastMarchingImageFilter::NodeContainer
typename LevelSetType::NodeContainer NodeContainer
Definition: itkFastMarchingImageFilter.h:131
itk::FastMarchingImageFilter::NodeContainerPointer
typename LevelSetType::NodeContainerPointer NodeContainerPointer
Definition: itkFastMarchingImageFilter.h:132
itk::FastMarchingImageFilter::LevelSetImageType
typename LevelSetType::LevelSetImageType LevelSetImageType
Definition: itkFastMarchingImageFilter.h:126
itk::FastMarchingUpwindGradientImageFilter::m_TargetPoints
NodeContainerPointer m_TargetPoints
Definition: itkFastMarchingUpwindGradientImageFilter.h:234
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::FastMarchingImageFilter::SpeedImageConstPointer
typename SpeedImageType::ConstPointer SpeedImageConstPointer
Definition: itkFastMarchingImageFilter.h:169
itk::FastMarchingUpwindGradientImageFilter::AxisNodeType
typename Superclass::AxisNodeType AxisNodeType
Definition: itkFastMarchingUpwindGradientImageFilter.h:85
itk::FastMarchingUpwindGradientImageFilter::m_TargetOffset
double m_TargetOffset
Definition: itkFastMarchingUpwindGradientImageFilter.h:241
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:83