ITK  5.3.0
Insight Toolkit
itkFastMarchingUpwindGradientImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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  * https://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 {
60 template <typename TLevelSet, typename TSpeedImage = Image<float, TLevelSet::ImageDimension>>
61 class ITK_TEMPLATE_EXPORT FastMarchingUpwindGradientImageFilter : public FastMarchingImageFilter<TLevelSet, TSpeedImage>
62 {
63 public:
64  ITK_DISALLOW_COPY_AND_MOVE(FastMarchingUpwindGradientImageFilter);
65 
71 
73  itkNewMacro(Self);
74 
77 
79  using typename Superclass::LevelSetType;
80  using typename Superclass::SpeedImageType;
81  using typename Superclass::LevelSetImageType;
82  using typename Superclass::LevelSetPointer;
83  using typename Superclass::SpeedImageConstPointer;
84  using typename Superclass::LabelImageType;
85  using typename Superclass::PixelType;
86  using typename Superclass::AxisNodeType;
87  using typename Superclass::NodeType;
88  using typename Superclass::NodeContainer;
89  using typename Superclass::NodeContainerPointer;
90 
91  using typename Superclass::IndexType;
92  using typename Superclass::OutputSpacingType;
93  using typename Superclass::LevelSetIndexType;
94 
95  using PointType = typename Superclass::OutputPointType;
96 
98  static constexpr unsigned int SetDimension = Superclass::SetDimension;
99 
103  void
105  {
106  m_TargetPoints = points;
107  this->Modified();
108  }
112  NodeContainerPointer
114  {
115  return m_TargetPoints;
116  }
117 
119  NodeContainerPointer
121  {
122  return m_ReachedTargetPoints;
123  }
124 
127 
130 
133 
137  {
138  return m_GradientImage;
139  }
140 
143  itkSetMacro(GenerateGradientImage, bool);
144 
146  itkGetConstReferenceMacro(GenerateGradientImage, bool);
147  itkBooleanMacro(GenerateGradientImage);
153  itkSetMacro(TargetOffset, double);
154 
156  itkGetConstReferenceMacro(TargetOffset, double);
157 
161  itkSetMacro(TargetReachedMode, int);
162  itkGetConstReferenceMacro(TargetReachedMode, int);
163  void
165  {
166  this->SetTargetReachedMode(NoTargets);
167  m_NumberOfTargets = 0;
168  }
169  void
171  {
172  this->VerifyTargetReachedModeConditions();
175  this->SetTargetReachedMode(OneTarget);
176  m_NumberOfTargets = 1;
177  }
178  void
180  {
181  this->VerifyTargetReachedModeConditions(numberOfTargets);
182 
183  this->SetTargetReachedMode(SomeTargets);
184  m_NumberOfTargets = numberOfTargets;
185  }
186 
187  void
189  {
190  this->VerifyTargetReachedModeConditions();
191 
192  this->SetTargetReachedMode(AllTargets);
193  m_NumberOfTargets = m_TargetPoints->Size();
194  }
195 
197  itkGetConstReferenceMacro(NumberOfTargets, SizeValueType);
198 
203  itkGetConstReferenceMacro(TargetValue, double);
204 
205  enum
206  {
210  AllTargets
211  };
212 
213 #ifdef ITK_USE_CONCEPT_CHECKING
214  // Begin concept checking
215  itkConceptMacro(LevelSetDoubleDivisionOperatorsCheck,
217  itkConceptMacro(LevelSetDoubleDivisionAndAssignOperatorsCheck,
219  // End concept checking
220 #endif
221 
222 protected:
224  ~FastMarchingUpwindGradientImageFilter() override = default;
225  void
226  PrintSelf(std::ostream & os, Indent indent) const override;
227 
228  void
229  Initialize(LevelSetImageType *) override;
230 
231  void
232  GenerateData() override;
233 
234  void
235  UpdateNeighbors(const IndexType & index, const SpeedImageType *, LevelSetImageType *) override;
236 
237  virtual void
238  ComputeGradient(const IndexType & index,
239  const LevelSetImageType * output,
240  const LabelImageType * labelImage,
241  GradientImageType * gradientImage);
242 
246  bool
248  {
249  if (!m_TargetPoints || m_TargetPoints->Size() == 0)
250  {
251  return false;
252  }
253  else
254  {
255  return true;
256  }
257  }
265  void
266  VerifyTargetReachedModeConditions(unsigned int targetModeMinPoints = 1)
267  {
268  bool targetPointsExist = this->IsTargetPointsExistenceConditionSatisfied();
271  if (!targetPointsExist)
272  {
273  itkExceptionMacro(<< "No target point set. Cannot set the target reached mode.");
274  }
275  else
276  {
277  SizeValueType availableNumberOfTargets = m_TargetPoints->Size();
278  if (targetModeMinPoints > availableNumberOfTargets)
279  {
280  itkExceptionMacro(<< "Not enough target points: Available: " << availableNumberOfTargets
281  << "; Requested: " << targetModeMinPoints);
282  }
283  }
284  }
285 
286 
287 private:
290 
292 
294 
296 
298 
300 
302 };
303 } // namespace itk
304 
305 #ifndef ITK_MANUAL_INSTANTIATION
306 # include "itkFastMarchingUpwindGradientImageFilter.hxx"
307 #endif
308 
309 #endif
itk::FastMarchingUpwindGradientImageFilter::OneTarget
Definition: itkFastMarchingUpwindGradientImageFilter.h:208
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:92
itk::FastMarchingUpwindGradientImageFilter::m_TargetReachedMode
int m_TargetReachedMode
Definition: itkFastMarchingUpwindGradientImageFilter.h:297
itk::FastMarchingUpwindGradientImageFilter::SetTargetReachedModeToAllTargets
void SetTargetReachedModeToAllTargets()
Definition: itkFastMarchingUpwindGradientImageFilter.h:188
itk::FastMarchingUpwindGradientImageFilter::m_GenerateGradientImage
bool m_GenerateGradientImage
Definition: itkFastMarchingUpwindGradientImageFilter.h:293
itk::FastMarchingUpwindGradientImageFilter::m_GradientImage
GradientImagePointer m_GradientImage
Definition: itkFastMarchingUpwindGradientImageFilter.h:291
itk::FastMarchingUpwindGradientImageFilter::NoTargets
Definition: itkFastMarchingUpwindGradientImageFilter.h:207
itkImage.h
itk::FastMarchingUpwindGradientImageFilter::m_ReachedTargetPoints
NodeContainerPointer m_ReachedTargetPoints
Definition: itkFastMarchingUpwindGradientImageFilter.h:289
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::FastMarchingUpwindGradientImageFilter::GetReachedTargetPoints
NodeContainerPointer GetReachedTargetPoints()
Definition: itkFastMarchingUpwindGradientImageFilter.h:120
itk::FastMarchingUpwindGradientImageFilter
Generates the upwind gradient field of fast marching arrival times.
Definition: itkFastMarchingUpwindGradientImageFilter.h:61
itk::Concept::DivisionAndAssignOperators
Definition: itkConceptChecking.h:502
itk::FastMarchingUpwindGradientImageFilter::GetTargetPoints
NodeContainerPointer GetTargetPoints()
Definition: itkFastMarchingUpwindGradientImageFilter.h:113
itk::FastMarchingUpwindGradientImageFilter::m_NumberOfTargets
SizeValueType m_NumberOfTargets
Definition: itkFastMarchingUpwindGradientImageFilter.h:301
itk::FastMarchingUpwindGradientImageFilter::SetTargetReachedModeToNoTargets
void SetTargetReachedModeToNoTargets()
Definition: itkFastMarchingUpwindGradientImageFilter.h:164
itkFastMarchingImageFilter.h
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::FastMarchingUpwindGradientImageFilter::GetGradientImage
GradientImagePointer GetGradientImage() const
Definition: itkFastMarchingUpwindGradientImageFilter.h:136
itk::FastMarchingUpwindGradientImageFilter::SetTargetPoints
void SetTargetPoints(NodeContainer *points)
Definition: itkFastMarchingUpwindGradientImageFilter.h:104
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itk::FastMarchingImageFilter
Solve an Eikonal equation using Fast Marching.
Definition: itkFastMarchingImageFilter.h:136
itk::FastMarchingUpwindGradientImageFilter::SetTargetReachedModeToSomeTargets
void SetTargetReachedModeToSomeTargets(SizeValueType numberOfTargets)
Definition: itkFastMarchingUpwindGradientImageFilter.h:179
itk::FastMarchingUpwindGradientImageFilter::m_TargetValue
double m_TargetValue
Definition: itkFastMarchingUpwindGradientImageFilter.h:299
itk::FastMarchingUpwindGradientImageFilter::VerifyTargetReachedModeConditions
void VerifyTargetReachedModeConditions(unsigned int targetModeMinPoints=1)
Definition: itkFastMarchingUpwindGradientImageFilter.h:266
itk::Concept::DivisionOperators
Definition: itkConceptChecking.h:473
itk::FastMarchingUpwindGradientImageFilter::SomeTargets
Definition: itkFastMarchingUpwindGradientImageFilter.h:209
itk::FastMarchingUpwindGradientImageFilter::PointType
typename Superclass::OutputPointType PointType
Definition: itkFastMarchingUpwindGradientImageFilter.h:95
itk::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition: itkCovariantVector.h:70
itk::FastMarchingUpwindGradientImageFilter::SetTargetReachedModeToOneTarget
void SetTargetReachedModeToOneTarget()
Definition: itkFastMarchingUpwindGradientImageFilter.h:170
itk::FastMarchingUpwindGradientImageFilter::GradientImagePointer
typename GradientImageType::Pointer GradientImagePointer
Definition: itkFastMarchingUpwindGradientImageFilter.h:132
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::FastMarchingImageFilter::NodeContainer
typename LevelSetType::NodeContainer NodeContainer
Definition: itkFastMarchingImageFilter.h:160
itk::FastMarchingImageFilter::NodeContainerPointer
typename LevelSetType::NodeContainerPointer NodeContainerPointer
Definition: itkFastMarchingImageFilter.h:161
itk::FastMarchingUpwindGradientImageFilter::m_TargetPoints
NodeContainerPointer m_TargetPoints
Definition: itkFastMarchingUpwindGradientImageFilter.h:288
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::FastMarchingUpwindGradientImageFilter::m_TargetOffset
double m_TargetOffset
Definition: itkFastMarchingUpwindGradientImageFilter.h:295
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:83
itk::FastMarchingUpwindGradientImageFilter::IsTargetPointsExistenceConditionSatisfied
bool IsTargetPointsExistenceConditionSatisfied()
Definition: itkFastMarchingUpwindGradientImageFilter.h:247