ITK  6.0.0
Insight Toolkit
itkPathToImageFilter.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 itkPathToImageFilter_h
19#define itkPathToImageFilter_h
20
21#include "itkImageSource.h"
22#include "itkConceptChecking.h"
23
24namespace itk
25{
36template <typename TInputPath, typename TOutputImage>
37class ITK_TEMPLATE_EXPORT PathToImageFilter : public ImageSource<TOutputImage>
38{
39public:
40 ITK_DISALLOW_COPY_AND_MOVE(PathToImageFilter);
41
47
49 itkNewMacro(Self);
50
52 itkOverrideGetNameOfClassMacro(PathToImageFilter);
53
55 using typename Superclass::OutputImageRegionType;
56 using InputPathType = TInputPath;
59 using OutputImageType = TOutputImage;
62 using ValueType = typename OutputImageType::ValueType;
63
65 static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
66
68 using Superclass::SetInput;
69 virtual void
70 SetInput(const InputPathType * input);
71
72 virtual void
73 SetInput(unsigned int, const TInputPath * path);
74
75 const InputPathType *
77
78 const InputPathType *
79 GetInput(unsigned int idx);
80
85 virtual void
86 SetSpacing(const double * spacing);
87
88 virtual void
89 SetSpacing(const float * spacing);
90
91 virtual const double *
92 GetSpacing() const;
93
96 itkSetMacro(PathValue, ValueType);
97 itkGetConstMacro(PathValue, ValueType);
98 itkSetMacro(BackgroundValue, ValueType);
99 itkGetConstMacro(BackgroundValue, ValueType);
106 virtual void
107 SetOrigin(const double * origin);
108
109 virtual void
110 SetOrigin(const float * origin);
111
112 virtual const double *
113 GetOrigin() const;
114
116 itkSetMacro(Size, SizeType);
117 itkGetConstMacro(Size, SizeType);
120protected:
122 ~PathToImageFilter() override = default;
123
124 void
126 {} // do nothing
127 void
128 GenerateData() override;
129
130 SizeType m_Size{};
131 double m_Spacing[OutputImageDimension]{};
132 double m_Origin[OutputImageDimension]{};
133 ValueType m_PathValue{};
134 ValueType m_BackgroundValue{};
135
136 void
137 PrintSelf(std::ostream & os, Indent indent) const override;
138};
139} // end namespace itk
140
141#ifndef ITK_MANUAL_INSTANTIATION
142# include "itkPathToImageFilter.hxx"
143#endif
144
145#endif
Base class for all process objects that output image data.
TOutputImage OutputImageType
typename OutputImageType::Pointer OutputImagePointer
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Base class for filters that take a Path as input and produce an image as output. Base class for filte...
void PrintSelf(std::ostream &os, Indent indent) const override
virtual void SetSpacing(const float *spacing)
typename InputPathType::ConstPointer InputPathConstPointer
void GenerateOutputInformation() override
const InputPathType * GetInput(unsigned int idx)
virtual void SetSpacing(const double *spacing)
typename InputPathType::Pointer InputPathPointer
virtual const double * GetSpacing() const
typename OutputImageType::ValueType ValueType
virtual void SetOrigin(const double *origin)
void GenerateData() override
~PathToImageFilter() override=default
virtual void SetInput(unsigned int, const TInputPath *path)
virtual const double * GetOrigin() const
const InputPathType * GetInput()
virtual void SetInput(const InputPathType *input)
virtual void SetOrigin(const float *origin)
typename OutputImageType::SizeType SizeType
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:70