ITK  6.0.0
Insight Toolkit
itkFastMarchingImageFilter.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 itkFastMarchingImageFilter_h
19#define itkFastMarchingImageFilter_h
20
23#include "itkLevelSet.h"
24#include "itkMath.h"
25#include "ITKFastMarchingExport.h"
26
27#include <functional>
28#include <queue>
29
30namespace itk
31{
37{
38public:
46 enum class Label : uint8_t
47 {
48 FarPoint = 0,
49 AlivePoint,
50 TrialPoint,
51 InitialTrialPoint,
52 OutsidePoint
53 };
54};
55// Define how to print enumeration
56extern ITKFastMarching_EXPORT std::ostream &
57 operator<<(std::ostream & out, const FastMarchingImageFilterEnums::Label value);
58
134template <typename TLevelSet, typename TSpeedImage = Image<float, TLevelSet::ImageDimension>>
135class ITK_TEMPLATE_EXPORT FastMarchingImageFilter : public ImageToImageFilter<TSpeedImage, TLevelSet>
136{
137public:
138 ITK_DISALLOW_COPY_AND_MOVE(FastMarchingImageFilter);
139
145
147 itkNewMacro(Self);
148
150 itkOverrideGetNameOfClassMacro(FastMarchingImageFilter);
151
163 using OutputSpacingType = typename LevelSetImageType::SpacingType;
166
167 class AxisNodeType : public NodeType
168 {
169 public:
170 AxisNodeType() = default;
171 int
172 GetAxis() const
173 {
174 return m_Axis;
175 }
176 void
177 SetAxis(int axis)
178 {
179 m_Axis = axis;
180 }
182 operator=(const NodeType & node)
183 {
184 this->NodeType::operator=(node);
185 return *this;
186 }
187
188 private:
189 int m_Axis{ 0 };
190 };
191
193 using SpeedImageType = TSpeedImage;
194
198
200 static constexpr unsigned int SetDimension = LevelSetType::SetDimension;
201 static constexpr unsigned int SpeedImageDimension = SpeedImageType::ImageDimension;
202
205
207#if !defined(ITK_LEGACY_REMOVE)
209 static constexpr LabelEnum FarPoint = LabelEnum::FarPoint;
210 static constexpr LabelEnum AlivePoint = LabelEnum::AlivePoint;
211 static constexpr LabelEnum TrialPoint = LabelEnum::TrialPoint;
212 static constexpr LabelEnum InitialTrialPoint = LabelEnum::InitialTrialPoint;
213 static constexpr LabelEnum OutsidePoint = LabelEnum::OutsidePoint;
214#endif
215
218
221
222 template <typename TPixel>
223 void
225 {
226 using InternalImageType = Image<TPixel, SetDimension>;
227 using InternalRegionIterator = ImageRegionConstIteratorWithIndex<InternalImageType>;
228 InternalRegionIterator b_it(iImage, iImage->GetLargestPossibleRegion());
229 b_it.GoToBegin();
230
231 TPixel zero_value{};
232 typename NodeContainer::ElementIdentifier NumberOfPoints = 0;
233
234 NodeType node;
235 node.SetValue(0.);
236
237 while (!b_it.IsAtEnd())
238 {
239 if (Math::ExactlyEquals(b_it.Get(), zero_value))
240 {
241 if (NumberOfPoints == 0)
242 {
243 m_OutsidePoints = NodeContainer::New();
244 }
245 node.SetIndex(b_it.GetIndex());
246 m_OutsidePoints->InsertElement(NumberOfPoints++, node);
247 }
248 ++b_it;
249 }
250 this->Modified();
251 }
252
254 void
256 {
257 m_OutsidePoints = points;
258 this->Modified();
259 }
264 void
266 {
267 m_AlivePoints = points;
268 this->Modified();
269 }
273 NodeContainerPointer
275 {
276 return m_AlivePoints;
277 }
278
281 void
283 {
284 m_TrialPoints = points;
285 this->Modified();
286 }
290 NodeContainerPointer
292 {
293 return m_TrialPoints;
294 }
295
297 LabelImagePointer
299 {
300 return m_LabelImage;
301 }
302
306 void
307 SetSpeedConstant(double value)
308 {
309 m_SpeedConstant = value;
310 m_InverseSpeed = -1.0 * itk::Math::sqr(1.0 / m_SpeedConstant);
311 this->Modified();
312 }
316 itkGetConstReferenceMacro(SpeedConstant, double);
317
322 itkSetMacro(NormalizationFactor, double);
323 itkGetConstMacro(NormalizationFactor, double);
329 itkSetMacro(StoppingValue, double);
330
332 itkGetConstReferenceMacro(StoppingValue, double);
333
338 itkSetMacro(CollectPoints, bool);
339
341 itkGetConstReferenceMacro(CollectPoints, bool);
342 itkBooleanMacro(CollectPoints);
351 {
352 return m_ProcessedPoints;
353 }
354
361 virtual void
363 {
364 m_OutputRegion = size;
365 }
366 virtual OutputSizeType
368 {
369 return m_OutputRegion.GetSize();
370 }
371 itkSetMacro(OutputRegion, OutputRegionType);
372 itkGetConstReferenceMacro(OutputRegion, OutputRegionType);
373 itkSetMacro(OutputSpacing, OutputSpacingType);
374 itkGetConstReferenceMacro(OutputSpacing, OutputSpacingType);
375 itkSetMacro(OutputDirection, OutputDirectionType);
376 itkGetConstReferenceMacro(OutputDirection, OutputDirectionType);
377 itkSetMacro(OutputOrigin, OutputPointType);
378 itkGetConstReferenceMacro(OutputOrigin, OutputPointType);
379 itkSetMacro(OverrideOutputInformation, bool);
380 itkGetConstReferenceMacro(OverrideOutputInformation, bool);
381 itkBooleanMacro(OverrideOutputInformation);
384#ifdef ITK_USE_CONCEPT_CHECKING
385 // Begin concept checking
388 itkConceptMacro(DoubleConvertibleToLevelSetCheck, (Concept::Convertible<double, PixelType>));
389 itkConceptMacro(LevelSetOStreamWritableCheck, (Concept::OStreamWritable<PixelType>));
390 // End concept checking
391#endif
392
393protected:
395 ~FastMarchingImageFilter() override = default;
396 void
397 PrintSelf(std::ostream & os, Indent indent) const override;
398
399 virtual void
401
402 virtual void
404
405 virtual double
407
408 const AxisNodeType &
409 GetNodeUsedInCalculation(unsigned int idx) const
410 {
411 return m_NodesUsed[idx];
412 }
413
414 void
415 GenerateData() override;
416
418 void
420
421 void
423
428 itkGetConstReferenceMacro(LargeValue, PixelType);
429
430 OutputRegionType m_BufferedRegion{};
432 LevelSetIndexType m_StartIndex{};
433 LevelSetIndexType m_LastIndex{};
434
435 itkGetConstReferenceMacro(StartIndex, LevelSetIndexType);
436 itkGetConstReferenceMacro(LastIndex, LevelSetIndexType);
437
438private:
439 NodeContainerPointer m_AlivePoints{};
440 NodeContainerPointer m_TrialPoints{};
441 NodeContainerPointer m_OutsidePoints{};
442
443 LabelImagePointer m_LabelImage{};
444
445 double m_SpeedConstant{};
446 double m_InverseSpeed{};
447 double m_StoppingValue{};
448
449 bool m_CollectPoints{};
450 NodeContainerPointer m_ProcessedPoints{};
451
452 OutputRegionType m_OutputRegion{};
453 OutputPointType m_OutputOrigin{};
454 OutputSpacingType m_OutputSpacing{};
455 OutputDirectionType m_OutputDirection{};
456 bool m_OverrideOutputInformation{};
457
458 typename LevelSetImageType::PixelType m_LargeValue{};
459 AxisNodeType m_NodesUsed[SetDimension]{};
460
464 using HeapContainer = std::vector<AxisNodeType>;
465 using NodeComparer = std::greater<AxisNodeType>;
466 using HeapType = std::priority_queue<AxisNodeType, HeapContainer, NodeComparer>;
467
468 HeapType m_TrialHeap{};
469
470 double m_NormalizationFactor{};
471};
472} // namespace itk
473
474#ifndef ITK_MANUAL_INSTANTIATION
475# include "itkFastMarchingImageFilter.hxx"
476#endif
477
478#endif
Base class for all data objects in ITK.
Contains all enum classes used by the FastMarchingImageFilter class.
AxisNodeType & operator=(const NodeType &node)
Solve an Eikonal equation using Fast Marching.
typename LevelSetType::NodeContainer NodeContainer
void GenerateOutputInformation() override
typename SpeedImageType::ConstPointer SpeedImageConstPointer
typename LevelSetImageType::SpacingType OutputSpacingType
const AxisNodeType & GetNodeUsedInCalculation(unsigned int idx) const
typename LevelSetType::NodeType NodeType
virtual OutputSizeType GetOutputSize() const
typename LevelSetImageType::SizeType OutputSizeType
virtual void Initialize(LevelSetImageType *)
typename LevelSetType::LevelSetPointer LevelSetPointer
std::vector< AxisNodeType > HeapContainer
typename LevelSetType::NodeContainerPointer NodeContainerPointer
void PrintSelf(std::ostream &os, Indent indent) const override
~FastMarchingImageFilter() override=default
void SetBinaryMask(Image< TPixel, SetDimension > *iImage)
void GenerateData() override
LabelImagePointer GetLabelImage() const
typename LevelSetImageType::RegionType OutputRegionType
typename LevelSetImageType::PointType OutputPointType
typename LevelSetImageType::IndexType LevelSetIndexType
std::priority_queue< AxisNodeType, HeapContainer, NodeComparer > HeapType
typename SpeedImageType::Pointer SpeedImagePointer
void SetAlivePoints(NodeContainer *points)
void SetOutsidePoints(NodeContainer *points)
typename LevelSetImageType::DirectionType OutputDirectionType
typename LevelSetType::LevelSetImageType LevelSetImageType
typename LevelSetType::PixelType PixelType
NodeContainerPointer GetProcessedPoints() const
void EnlargeOutputRequestedRegion(DataObject *output) override
typename LabelImageType::Pointer LabelImagePointer
void SetTrialPoints(NodeContainer *points)
virtual void UpdateNeighbors(const IndexType &index, const SpeedImageType *, LevelSetImageType *)
typename NodeType::IndexType NodeIndexType
virtual double UpdateValue(const IndexType &index, const SpeedImageType *, LevelSetImageType *)
virtual void SetOutputSize(const OutputSizeType &size)
std::greater< AxisNodeType > NodeComparer
virtual const RegionType & GetLargestPossibleRegion() const
Definition: itkImageBase.h:279
A multi-dimensional iterator templated over image type that walks an image region and is specialized ...
Base class for filters that take an image as input and produce an image as output.
Templated n-dimensional image class.
Definition: itkImage.h:89
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Represent a node in a level set.
Level set type information.
Definition: itkLevelSet.h:41
typename TLevelSet::PixelType PixelType
Definition: itkLevelSet.h:55
typename NodeContainer::Pointer NodeContainerPointer
Definition: itkLevelSet.h:64
typename TLevelSet::Pointer LevelSetPointer
Definition: itkLevelSet.h:51
Light weight base class for most itk classes.
Define a front-end to the STL "vector" container that conforms to the IndexedContainerInterface.
SmartPointer< const Self > ConstPointer
static Pointer New()
SmartPointer< Self > Pointer
#define itkConceptMacro(name, concept)
bool ExactlyEquals(const TInput1 &x1, const TInput2 &x2)
Return the result of an exact comparison between two scalar values of potentially different types.
Definition: itkMath.h:722
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)
Represent a n-dimensional index in a n-dimensional image.
Definition: itkIndex.h:69