ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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
48
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;
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;
117
118 using typename Superclass::IndexType;
119 using typename Superclass::OutputSpacingType;
120 using typename Superclass::LevelSetIndexType;
121
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 }
147
148
150 NodeContainerPointer
152 {
153 return m_TargetPoints;
154 }
155
157 NodeContainerPointer
162
165
168
171
175 {
176 return m_GradientImage;
177 }
178
181 itkSetMacro(GenerateGradientImage, bool);
182
184 itkGetConstReferenceMacro(GenerateGradientImage, bool);
185 itkBooleanMacro(GenerateGradientImage);
187
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);
206 }
207 void
209 {
210 this->SetTargetReachedMode(TargetConditionEnum::OneTarget);
212 }
213 void
215 {
216 this->SetTargetReachedMode(TargetConditionEnum::SomeTargets);
217 m_NumberOfTargets = numberOfTargets;
218 }
219
220
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 itkConceptMacro(LevelSetDoubleDivisionOperatorsCheck,
239 itkConceptMacro(LevelSetDoubleDivisionAndAssignOperatorsCheck,
241
242protected:
245 void
246 PrintSelf(std::ostream & os, Indent indent) const override;
247
248 void
249 VerifyPreconditions() const override;
250
251 void
253
254 void
255 GenerateData() override;
256
257 void
258 UpdateNeighbors(const IndexType & index, const SpeedImageType *, LevelSetImageType *) override;
259
260 virtual void
262 const LevelSetImageType * output,
263 const LabelImageType * labelImage,
264 GradientImageType * gradientImage);
265
269 bool
271 {
272 if (!m_TargetPoints || m_TargetPoints->Size() == 0)
273 {
274 return false;
275 }
277
278 return true;
279 }
280
288 void
289 VerifyTargetReachedModeConditions(unsigned int targetModeMinPoints = 1) const
290 {
291 const bool targetPointsExist = this->IsTargetPointsExistenceConditionSatisfied();
293
294 if (!targetPointsExist)
295 {
296 itkExceptionMacro("No target point set. Cannot set the target reached mode.");
297 }
298 else
299 {
300 const SizeValueType availableNumberOfTargets = m_TargetPoints->Size();
301 if (targetModeMinPoints > availableNumberOfTargets)
302 {
303 itkExceptionMacro("Not enough target points: Available: " << availableNumberOfTargets
304 << "; Requested: " << targetModeMinPoints);
305 }
306 }
307 }
308
309
310private:
313
315
317
319
321
323
325};
326} // namespace itk
327
328#ifndef ITK_MANUAL_INSTANTIATION
329# include "itkFastMarchingUpwindGradientImageFilter.hxx"
330#endif
331
332#endif
A templated class holding a n-Dimensional covariant vector.
typename LevelSetType::NodeContainer NodeContainer
typename SpeedImageType::ConstPointer SpeedImageConstPointer
typename LevelSetImageType::SpacingType OutputSpacingType
typename LevelSetType::NodeType NodeType
typename LevelSetType::LevelSetPointer LevelSetPointer
Image< LabelEnum, Self::SetDimension > LabelImageType
typename LevelSetType::NodeContainerPointer NodeContainerPointer
static constexpr unsigned int SetDimension
typename LevelSetImageType::IndexType LevelSetIndexType
typename LevelSetType::LevelSetImageType LevelSetImageType
LevelSetTypeDefault< TLevelSet > LevelSetType
Index< Self::SetDimension > IndexType
typename LevelSetType::PixelType PixelType
enums for itk::FastMarchingUpwindGradientImageFilter
Image< GradientPixelType, Self::SetDimension > GradientImageType
void PrintSelf(std::ostream &os, Indent indent) const override
void VerifyPreconditions() const override
Verifies that the process object has been configured correctly, that all required inputs are set,...
CovariantVector< PixelType, Self::SetDimension > GradientPixelType
virtual void ComputeGradient(const IndexType &index, const LevelSetImageType *output, const LabelImageType *labelImage, GradientImageType *gradientImage)
~FastMarchingUpwindGradientImageFilter() override=default
void Initialize(LevelSetImageType *) override
FastMarchingUpwindGradientImageFilterEnums::TargetCondition TargetConditionEnum
void VerifyTargetReachedModeConditions(unsigned int targetModeMinPoints=1) const
void UpdateNeighbors(const IndexType &index, const SpeedImageType *, LevelSetImageType *) override
FastMarchingImageFilter< TLevelSet, TSpeedImage > Superclass
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
virtual void Modified() const
Implements transparent reference counting.
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition itkIntTypes.h:86
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)