ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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{
55// Define how to print enumeration
56extern ITKFastMarching_EXPORT std::ostream &
57 operator<<(std::ostream & out, const FastMarchingImageFilterEnums::Label value);
58
132template <typename TLevelSet, typename TSpeedImage = Image<float, TLevelSet::ImageDimension>>
133class ITK_TEMPLATE_EXPORT FastMarchingImageFilter : public ImageToImageFilter<TSpeedImage, TLevelSet>
134{
135public:
136 ITK_DISALLOW_COPY_AND_MOVE(FastMarchingImageFilter);
137
143
145 itkNewMacro(Self);
146
148 itkOverrideGetNameOfClassMacro(FastMarchingImageFilter);
149
156 using NodeIndexType = typename NodeType::IndexType;
159 using OutputSizeType = typename LevelSetImageType::SizeType;
160 using OutputRegionType = typename LevelSetImageType::RegionType;
161 using OutputSpacingType = typename LevelSetImageType::SpacingType;
162 using OutputDirectionType = typename LevelSetImageType::DirectionType;
163 using OutputPointType = typename LevelSetImageType::PointType;
164
165 class AxisNodeType : public NodeType
166 {
167 public:
168 AxisNodeType() = default;
169 int
170 GetAxis() const
171 {
172 return m_Axis;
173 }
174 void
175 SetAxis(int axis)
176 {
177 m_Axis = axis;
178 }
180 operator=(const NodeType & node)
181 {
182 this->NodeType::operator=(node);
183 return *this;
184 }
185
186 private:
187 int m_Axis{ 0 };
188 };
189
191 using SpeedImageType = TSpeedImage;
192
194 using SpeedImagePointer = typename SpeedImageType::Pointer;
195 using SpeedImageConstPointer = typename SpeedImageType::ConstPointer;
196
198 static constexpr unsigned int SetDimension = LevelSetType::SetDimension;
199 static constexpr unsigned int SpeedImageDimension = SpeedImageType::ImageDimension;
200
203
205#if !defined(ITK_LEGACY_REMOVE)
207 static constexpr LabelEnum FarPoint = LabelEnum::FarPoint;
208 static constexpr LabelEnum AlivePoint = LabelEnum::AlivePoint;
209 static constexpr LabelEnum TrialPoint = LabelEnum::TrialPoint;
210 static constexpr LabelEnum InitialTrialPoint = LabelEnum::InitialTrialPoint;
211 static constexpr LabelEnum OutsidePoint = LabelEnum::OutsidePoint;
212#endif
213
216
219
220 template <typename TPixel>
221 void
223 {
224 using InternalImageType = Image<TPixel, SetDimension>;
225 using InternalRegionIterator = ImageRegionConstIteratorWithIndex<InternalImageType>;
226 InternalRegionIterator b_it(iImage, iImage->GetLargestPossibleRegion());
227 b_it.GoToBegin();
228
229 TPixel zero_value{};
230 typename NodeContainer::ElementIdentifier NumberOfPoints = 0;
231
232 NodeType node;
233 node.SetValue(0.);
234
235 while (!b_it.IsAtEnd())
236 {
237 if (Math::ExactlyEquals(b_it.Get(), zero_value))
238 {
239 if (NumberOfPoints == 0)
240 {
241 m_OutsidePoints = NodeContainer::New();
242 }
243 node.SetIndex(b_it.GetIndex());
244 m_OutsidePoints->InsertElement(NumberOfPoints++, node);
245 }
246 ++b_it;
247 }
248 this->Modified();
249 }
250
252 void
254 {
255 m_OutsidePoints = points;
256 this->Modified();
257 }
258
261 void
263 {
264 m_AlivePoints = points;
265 this->Modified();
266 }
267
269 NodeContainerPointer
271 {
272 return m_AlivePoints;
273 }
274
277 void
279 {
280 m_TrialPoints = points;
281 this->Modified();
282 }
283
285 NodeContainerPointer
287 {
288 return m_TrialPoints;
289 }
290
292 LabelImagePointer
294 {
295 return m_LabelImage;
296 }
297
301 void
302 SetSpeedConstant(double value)
303 {
304 m_SpeedConstant = value;
305 m_InverseSpeed = -1.0 * itk::Math::sqr(1.0 / m_SpeedConstant);
306 this->Modified();
307 }
308
310 itkGetConstReferenceMacro(SpeedConstant, double);
311
317 itkSetMacro(NormalizationFactor, double);
318 itkGetConstMacro(NormalizationFactor, double);
320
324 itkSetMacro(StoppingValue, double);
325
327 itkGetConstReferenceMacro(StoppingValue, double);
328
333 itkSetMacro(CollectPoints, bool);
334
337 itkGetConstReferenceMacro(CollectPoints, bool);
338 itkBooleanMacro(CollectPoints);
340
347 {
348 return m_ProcessedPoints;
349 }
350
358 virtual void
360 {
361 m_OutputRegion = size;
362 }
363 virtual OutputSizeType
365 {
366 return m_OutputRegion.GetSize();
367 }
368 itkSetMacro(OutputRegion, OutputRegionType);
369 itkGetConstReferenceMacro(OutputRegion, OutputRegionType);
370 itkSetMacro(OutputSpacing, OutputSpacingType);
371 itkGetConstReferenceMacro(OutputSpacing, OutputSpacingType);
372 itkSetMacro(OutputDirection, OutputDirectionType);
373 itkGetConstReferenceMacro(OutputDirection, OutputDirectionType);
374 itkSetMacro(OutputOrigin, OutputPointType);
375 itkGetConstReferenceMacro(OutputOrigin, OutputPointType);
376 itkSetMacro(OverrideOutputInformation, bool);
377 itkGetConstReferenceMacro(OverrideOutputInformation, bool);
378 itkBooleanMacro(OverrideOutputInformation);
380
383 itkConceptMacro(DoubleConvertibleToLevelSetCheck, (Concept::Convertible<double, PixelType>));
384 itkConceptMacro(LevelSetOStreamWritableCheck, (Concept::OStreamWritable<PixelType>));
385
386protected:
388 ~FastMarchingImageFilter() override = default;
389 void
390 PrintSelf(std::ostream & os, Indent indent) const override;
391
392 virtual void
394
395 virtual void
397
398 virtual double
400
401 const AxisNodeType &
402 GetNodeUsedInCalculation(unsigned int idx) const
403 {
404 return m_NodesUsed[idx];
405 }
406
407 void
408 GenerateData() override;
409
411 void
413
414 void
416
421 itkGetConstReferenceMacro(LargeValue, PixelType);
422
424 using LevelSetIndexType = typename LevelSetImageType::IndexType;
427
428 itkGetConstReferenceMacro(StartIndex, LevelSetIndexType);
429 itkGetConstReferenceMacro(LastIndex, LevelSetIndexType);
430
431private:
435
437
441
444
450
451 typename LevelSetImageType::PixelType m_LargeValue{};
452 AxisNodeType m_NodesUsed[SetDimension]{};
453
457 using HeapContainer = std::vector<AxisNodeType>;
458 using NodeComparer = std::greater<AxisNodeType>;
459 using HeapType = std::priority_queue<AxisNodeType, HeapContainer, NodeComparer>;
460
462
464};
465} // namespace itk
466
467#ifndef ITK_MANUAL_INSTANTIATION
468# include "itkFastMarchingImageFilter.hxx"
469#endif
470
471#endif
Base class for all data objects in ITK.
Contains all enum classes used by the FastMarchingImageFilter class.
AxisNodeType & operator=(const NodeType &node)
void GenerateOutputInformation() override
const AxisNodeType & GetNodeUsedInCalculation(unsigned int idx) const
virtual OutputSizeType GetOutputSize() const
virtual void Initialize(LevelSetImageType *)
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
std::priority_queue< AxisNodeType, HeapContainer, NodeComparer > HeapType
void SetAlivePoints(NodeContainer *points)
void SetOutsidePoints(NodeContainer *points)
typename LevelSetImageType::DirectionType OutputDirectionType
ImageToImageFilter< TSpeedImage, TLevelSet > Superclass
NodeContainerPointer GetProcessedPoints() const
void EnlargeOutputRequestedRegion(DataObject *output) override
void SetTrialPoints(NodeContainer *points)
virtual void UpdateNeighbors(const IndexType &index, const SpeedImageType *, LevelSetImageType *)
virtual double UpdateValue(const IndexType &index, const SpeedImageType *, LevelSetImageType *)
virtual void SetOutputSize(const OutputSizeType &size)
virtual const RegionType & GetLargestPossibleRegion() const
A multi-dimensional iterator templated over image type that walks an image region and is specialized ...
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
Level set type information.
Definition itkLevelSet.h:41
LevelSetNode< PixelType, Self::SetDimension > NodeType
Definition itkLevelSet.h:58
typename TLevelSet::PixelType PixelType
Definition itkLevelSet.h:55
typename NodeContainer::Pointer NodeContainerPointer
Definition itkLevelSet.h:64
static constexpr unsigned int SetDimension
Definition itkLevelSet.h:48
VectorContainer< unsigned int, NodeType > NodeContainer
Definition itkLevelSet.h:61
typename TLevelSet::Pointer LevelSetPointer
Definition itkLevelSet.h:51
virtual void Modified() const
Implements transparent reference counting.
#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:719
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