ITK  6.0.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 {
26 
34 {
35 public:
40  enum class TargetCondition : int
41  {
42  NoTargets,
43  OneTarget,
46  };
47 };
49 extern ITKFastMarching_EXPORT std::ostream &
51 
52 
87 template <typename TLevelSet, typename TSpeedImage = Image<float, TLevelSet::ImageDimension>>
88 class ITK_TEMPLATE_EXPORT FastMarchingUpwindGradientImageFilter : public FastMarchingImageFilter<TLevelSet, TSpeedImage>
89 {
90 public:
91  ITK_DISALLOW_COPY_AND_MOVE(FastMarchingUpwindGradientImageFilter);
92 
98 
100  itkNewMacro(Self);
101 
103  itkOverrideGetNameOfClassMacro(FastMarchingUpwindGradientImageFilter);
104 
106  using typename Superclass::LevelSetType;
107  using typename Superclass::SpeedImageType;
108  using typename Superclass::LevelSetImageType;
109  using typename Superclass::LevelSetPointer;
110  using typename Superclass::SpeedImageConstPointer;
111  using typename Superclass::LabelImageType;
112  using typename Superclass::PixelType;
113  using typename Superclass::AxisNodeType;
114  using typename Superclass::NodeType;
115  using typename Superclass::NodeContainer;
116  using typename Superclass::NodeContainerPointer;
117 
118  using typename Superclass::IndexType;
119  using typename Superclass::OutputSpacingType;
120  using typename Superclass::LevelSetIndexType;
121 
122  using PointType = typename Superclass::OutputPointType;
123 
125  static constexpr unsigned int SetDimension = Superclass::SetDimension;
126 
129 #if !defined(ITK_LEGACY_REMOVE)
130  // We need to expose the enum values at the class level
131  // for backwards compatibility
132  static constexpr TargetConditionEnum NoTargets = TargetConditionEnum::NoTargets;
133  static constexpr TargetConditionEnum OneTarget = TargetConditionEnum::OneTarget;
134  static constexpr TargetConditionEnum SomeTargets = TargetConditionEnum::SomeTargets;
135  static constexpr TargetConditionEnum AllTargets = TargetConditionEnum::AllTargets;
136 #endif
137 
141  void
143  {
144  m_TargetPoints = points;
145  this->Modified();
146  }
150  NodeContainerPointer
152  {
153  return m_TargetPoints;
154  }
155 
157  NodeContainerPointer
159  {
160  return m_ReachedTargetPoints;
161  }
162 
165 
168 
171 
175  {
176  return m_GradientImage;
177  }
178 
181  itkSetMacro(GenerateGradientImage, bool);
182 
184  itkGetConstReferenceMacro(GenerateGradientImage, bool);
185  itkBooleanMacro(GenerateGradientImage);
191  itkSetMacro(TargetOffset, double);
192 
194  itkGetConstReferenceMacro(TargetOffset, double);
195 
199  itkSetMacro(TargetReachedMode, TargetConditionEnum);
200  itkGetConstReferenceMacro(TargetReachedMode, TargetConditionEnum);
201  void
203  {
204  this->SetTargetReachedMode(TargetConditionEnum::NoTargets);
205  m_NumberOfTargets = 0;
206  }
207  void
209  {
210  this->SetTargetReachedMode(TargetConditionEnum::OneTarget);
211  m_NumberOfTargets = 1;
212  }
213  void
215  {
216  this->SetTargetReachedMode(TargetConditionEnum::SomeTargets);
217  m_NumberOfTargets = numberOfTargets;
218  }
221  void
223  {
224  this->SetTargetReachedMode(TargetConditionEnum::AllTargets);
225  // m_NumberOfTargets is not used for this case
226  }
227 
229  itkGetConstReferenceMacro(NumberOfTargets, SizeValueType);
230 
235  itkGetConstReferenceMacro(TargetValue, double);
236 
237 #ifdef ITK_USE_CONCEPT_CHECKING
238  // Begin concept checking
239  itkConceptMacro(LevelSetDoubleDivisionOperatorsCheck,
241  itkConceptMacro(LevelSetDoubleDivisionAndAssignOperatorsCheck,
243  // End concept checking
244 #endif
245 
246 protected:
248  ~FastMarchingUpwindGradientImageFilter() override = default;
249  void
250  PrintSelf(std::ostream & os, Indent indent) const override;
251 
252  virtual void
253  VerifyPreconditions() const override;
254 
255  void
256  Initialize(LevelSetImageType *) override;
257 
258  void
259  GenerateData() override;
260 
261  void
262  UpdateNeighbors(const IndexType & index, const SpeedImageType *, LevelSetImageType *) override;
263 
264  virtual void
265  ComputeGradient(const IndexType & index,
266  const LevelSetImageType * output,
267  const LabelImageType * labelImage,
268  GradientImageType * gradientImage);
269 
273  bool
275  {
276  if (!m_TargetPoints || m_TargetPoints->Size() == 0)
277  {
278  return false;
279  }
280  else
281  {
282  return true;
283  }
284  }
294  void
295  VerifyTargetReachedModeConditions(unsigned int targetModeMinPoints = 1) const
296  {
297  const bool targetPointsExist = this->IsTargetPointsExistenceConditionSatisfied();
300  if (!targetPointsExist)
301  {
302  itkExceptionMacro("No target point set. Cannot set the target reached mode.");
303  }
304  else
305  {
306  const SizeValueType availableNumberOfTargets = m_TargetPoints->Size();
307  if (targetModeMinPoints > availableNumberOfTargets)
308  {
309  itkExceptionMacro("Not enough target points: Available: " << availableNumberOfTargets
310  << "; Requested: " << targetModeMinPoints);
311  }
312  }
313  }
314 
315 
316 private:
317  NodeContainerPointer m_TargetPoints{};
318  NodeContainerPointer m_ReachedTargetPoints{};
319 
320  GradientImagePointer m_GradientImage{};
321 
322  bool m_GenerateGradientImage{};
323 
324  double m_TargetOffset{};
325 
326  TargetConditionEnum m_TargetReachedMode{};
327 
328  double m_TargetValue{};
329 
330  SizeValueType m_NumberOfTargets{};
331 };
332 } // namespace itk
333 
334 #ifndef ITK_MANUAL_INSTANTIATION
335 # include "itkFastMarchingUpwindGradientImageFilter.hxx"
336 #endif
337 
338 #endif
itk::FastMarchingUpwindGradientImageFilterEnums
enums for itk::FastMarchingUpwindGradientImageFilter
Definition: itkFastMarchingUpwindGradientImageFilter.h:33
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::FastMarchingUpwindGradientImageFilter::SetTargetReachedModeToAllTargets
void SetTargetReachedModeToAllTargets()
Definition: itkFastMarchingUpwindGradientImageFilter.h:222
itk::FastMarchingUpwindGradientImageFilterEnums::TargetCondition::OneTarget
itk::FastMarchingUpwindGradientImageFilterEnums::TargetCondition::SomeTargets
itk::FastMarchingUpwindGradientImageFilter::IsTargetPointsExistenceConditionSatisfied
bool IsTargetPointsExistenceConditionSatisfied() const
Definition: itkFastMarchingUpwindGradientImageFilter.h:274
itkImage.h
itk::operator<<
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::FastMarchingUpwindGradientImageFilter::GetReachedTargetPoints
NodeContainerPointer GetReachedTargetPoints()
Definition: itkFastMarchingUpwindGradientImageFilter.h:158
itk::FastMarchingUpwindGradientImageFilter
Generates the upwind gradient field of fast marching arrival times.
Definition: itkFastMarchingUpwindGradientImageFilter.h:88
itk::Concept::DivisionAndAssignOperators
Definition: itkConceptChecking.h:505
itk::FastMarchingUpwindGradientImageFilter::GetTargetPoints
NodeContainerPointer GetTargetPoints()
Definition: itkFastMarchingUpwindGradientImageFilter.h:151
itk::FastMarchingUpwindGradientImageFilter::SetTargetReachedModeToNoTargets
void SetTargetReachedModeToNoTargets()
Definition: itkFastMarchingUpwindGradientImageFilter.h:202
itkFastMarchingImageFilter.h
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::FastMarchingUpwindGradientImageFilter::GetGradientImage
GradientImagePointer GetGradientImage() const
Definition: itkFastMarchingUpwindGradientImageFilter.h:174
itk::FastMarchingUpwindGradientImageFilter::SetTargetPoints
void SetTargetPoints(NodeContainer *points)
Definition: itkFastMarchingUpwindGradientImageFilter.h:142
itk::FastMarchingUpwindGradientImageFilter::VerifyTargetReachedModeConditions
void VerifyTargetReachedModeConditions(unsigned int targetModeMinPoints=1) const
Definition: itkFastMarchingUpwindGradientImageFilter.h:295
itk::FastMarchingUpwindGradientImageFilterEnums::TargetCondition::NoTargets
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
TargetCondition
itk::FastMarchingImageFilter
Solve an Eikonal equation using Fast Marching.
Definition: itkFastMarchingImageFilter.h:135
itk::FastMarchingUpwindGradientImageFilter::SetTargetReachedModeToSomeTargets
void SetTargetReachedModeToSomeTargets(SizeValueType numberOfTargets)
Definition: itkFastMarchingUpwindGradientImageFilter.h:214
itk::Concept::DivisionOperators
Definition: itkConceptChecking.h:476
itk::FastMarchingUpwindGradientImageFilter::PointType
typename Superclass::OutputPointType PointType
Definition: itkFastMarchingUpwindGradientImageFilter.h:122
itk::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition: itkCovariantVector.h:70
itk::FastMarchingUpwindGradientImageFilterEnums::TargetCondition
TargetCondition
Definition: itkFastMarchingUpwindGradientImageFilter.h:40
itk::FastMarchingUpwindGradientImageFilter::SetTargetReachedModeToOneTarget
void SetTargetReachedModeToOneTarget()
Definition: itkFastMarchingUpwindGradientImageFilter.h:208
itk::FastMarchingUpwindGradientImageFilter::GradientImagePointer
typename GradientImageType::Pointer GradientImagePointer
Definition: itkFastMarchingUpwindGradientImageFilter.h:170
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::FastMarchingImageFilter::NodeContainer
typename LevelSetType::NodeContainer NodeContainer
Definition: itkFastMarchingImageFilter.h:159
itk::FastMarchingImageFilter::NodeContainerPointer
typename LevelSetType::NodeContainerPointer NodeContainerPointer
Definition: itkFastMarchingImageFilter.h:160
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::FastMarchingUpwindGradientImageFilterEnums::TargetCondition::AllTargets
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:86