ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkPathToImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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  * http://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 
24 namespace itk
25 {
35 template< typename TInputPath, typename TOutputImage >
36 class ITK_TEMPLATE_EXPORT PathToImageFilter:public ImageSource< TOutputImage >
37 {
38 public:
44 
46  itkNewMacro(Self);
47 
50 
52  typedef typename Superclass::OutputImageRegionType OutputImageRegionType;
53  typedef TInputPath InputPathType;
54  typedef typename InputPathType::Pointer InputPathPointer;
55  typedef typename InputPathType::ConstPointer InputPathConstPointer;
56  typedef TOutputImage OutputImageType;
57  typedef typename OutputImageType::Pointer OutputImagePointer;
59  typedef typename OutputImageType::ValueType ValueType;
60 
62  itkStaticConstMacro(OutputImageDimension, unsigned int,
63  TOutputImage::ImageDimension);
64 
66  using Superclass::SetInput;
67  virtual void SetInput(const InputPathType *path);
68 
69  virtual void SetInput(unsigned int, const TInputPath *path);
70 
71  const InputPathType * GetInput();
72 
73  const InputPathType * GetInput(unsigned int idx);
74 
79  virtual void SetSpacing(const double *spacing);
80 
81  virtual void SetSpacing(const float *spacing);
82 
83  virtual const double * GetSpacing() const;
84 
87  itkSetMacro(PathValue, ValueType);
88  itkGetConstMacro(PathValue, ValueType);
89  itkSetMacro(BackgroundValue, ValueType);
90  itkGetConstMacro(BackgroundValue, ValueType);
92 
97  virtual void SetOrigin(const double *origin);
98 
99  virtual void SetOrigin(const float *origin);
100 
101  virtual const double * GetOrigin() const;
102 
104  itkSetMacro(Size, SizeType);
105  itkGetConstMacro(Size, SizeType);
107 
108 protected:
110  ~PathToImageFilter() ITK_OVERRIDE;
111 
112  virtual void GenerateOutputInformation() ITK_OVERRIDE {} // do nothing
113  virtual void GenerateData() ITK_OVERRIDE;
114 
115  SizeType m_Size;
116  double m_Spacing[OutputImageDimension];
117  double m_Origin[OutputImageDimension];
118  ValueType m_PathValue;
119  ValueType m_BackgroundValue;
120 
121  virtual void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
122 
123 private:
124  ITK_DISALLOW_COPY_AND_ASSIGN(PathToImageFilter);
125 };
126 } // end namespace itk
127 
128 #ifndef ITK_MANUAL_INSTANTIATION
129 #include "itkPathToImageFilter.hxx"
130 #endif
131 
132 #endif
SmartPointer< Self > Pointer
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
Base class for all process objects that output image data.
OutputImageType::Pointer OutputImagePointer
Base class for filters that take a Path as input and produce an image as output. Base class for filte...
SmartPointer< const Self > ConstPointer
InputPathType::ConstPointer InputPathConstPointer
OutputImageType::SizeType SizeType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
InputPathType::Pointer InputPathPointer
Superclass::OutputImageRegionType OutputImageRegionType
OutputImageType::ValueType ValueType
ImageSource< TOutputImage > Superclass