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
24namespace itk
25{
26
34{
35public:
40 enum class TargetCondition : int
41 {
42 NoTargets,
43 OneTarget,
44 SomeTargets,
45 AllTargets
46 };
47};
49extern ITKFastMarching_EXPORT std::ostream &
51
52
87template <typename TLevelSet, typename TSpeedImage = Image<float, TLevelSet::ImageDimension>>
88class ITK_TEMPLATE_EXPORT FastMarchingUpwindGradientImageFilter : public FastMarchingImageFilter<TLevelSet, TSpeedImage>
89{
90public:
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
246protected:
249 void
250 PrintSelf(std::ostream & os, Indent indent) const override;
251
252 virtual void
253 VerifyPreconditions() const override;
254
255 void
257
258 void
259 GenerateData() override;
260
261 void
262 UpdateNeighbors(const IndexType & index, const SpeedImageType *, LevelSetImageType *) override;
263
264 virtual void
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 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 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
316private:
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
A templated class holding a n-Dimensional covariant vector.
Solve an Eikonal equation using Fast Marching.
enums for itk::FastMarchingUpwindGradientImageFilter
Generates the upwind gradient field of fast marching arrival times.
virtual void VerifyPreconditions() const override
Verifies that the process object has been configured correctly, that all required inputs are set,...
void PrintSelf(std::ostream &os, Indent indent) const override
virtual void ComputeGradient(const IndexType &index, const LevelSetImageType *output, const LabelImageType *labelImage, GradientImageType *gradientImage)
~FastMarchingUpwindGradientImageFilter() override=default
void Initialize(LevelSetImageType *) override
void VerifyTargetReachedModeConditions(unsigned int targetModeMinPoints=1) const
void UpdateNeighbors(const IndexType &index, const SpeedImageType *, LevelSetImageType *) override
Templated n-dimensional image class.
Definition: itkImage.h:89
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
SmartPointer< Self > Pointer
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
std::ostream & operator<<(std::ostream &os, const Array< TValue > &arr)
Definition: itkArray.h:216
unsigned long SizeValueType
Definition: itkIntTypes.h:86
Represent a n-dimensional index in a n-dimensional image.
Definition: itkIndex.h:69