ITK  4.8.0
Insight Segmentation and Registration Toolkit
itkInterpolateImageFilter.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 itkInterpolateImageFilter_h
19 #define itkInterpolateImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
23 
24 namespace itk
25 {
43 template< typename TInputImage, typename TOutputImage >
45  public ImageToImageFilter< TInputImage, TOutputImage >
46 {
47 public:
53 
55  itkNewMacro(Self);
56 
59 
66 
68  itkStaticConstMacro(ImageDimension, unsigned int,
69  TOutputImage::ImageDimension);
70  itkStaticConstMacro(IntermediateImageDimension, unsigned int,
71  TOutputImage::ImageDimension + 1);
73 
75  typedef typename TInputImage::PixelType InputPixelType;
78 
80  void SetInput1(const InputImageType *image)
81  { this->SetInput(image); }
83  { return this->GetInput(); }
85 
87  void SetInput2(const InputImageType *image);
88 
89  const InputImageType * GetInput2();
90 
93  itkSetClampMacro(Distance, double, 0.0, 1.0);
94  itkGetConstMacro(Distance, double);
96 
98  itkSetObjectMacro(Interpolator, InterpolatorType)
99  itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
100 
103  void BeforeThreadedGenerateData() ITK_OVERRIDE;
104 
106  void AfterThreadedGenerateData() ITK_OVERRIDE;
107 
108 #ifdef ITK_USE_CONCEPT_CHECKING
109  // Begin concept checking
110  itkConceptMacro( InputHasNumericTraitsCheck,
112  // End concept checking
113 #endif
114 
115 protected:
118  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
119 
121  void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
122  ThreadIdType threadId) ITK_OVERRIDE;
123 
124 private:
125  InterpolateImageFilter(const Self &); //purposely not implemented
126  void operator=(const Self &); //purposely not implemented
127 
129 
131 
132  double m_Distance;
133 };
134 } // end namespace itk
135 
136 #ifndef ITK_MANUAL_INSTANTIATION
137 #include "itkInterpolateImageFilter.hxx"
138 #endif
139 
140 #endif
Superclass::OutputImageType OutputImageType
void AfterThreadedGenerateData() override
Superclass::OutputImagePointer OutputImagePointer
static const unsigned int ImageDimension
InputImageType::Pointer InputImagePointer
static const unsigned int IntermediateImageDimension
void operator=(const Self &)
void BeforeThreadedGenerateData() override
const InputImageType * GetInput1()
SmartPointer< const Self > ConstPointer
Base class for all process objects that output image data.
Image< InputPixelType, itkGetStaticConstMacro(IntermediateImageDimension) > IntermediateImageType
virtual void SetInput(const InputImageType *image)
OutputImageType::Pointer OutputImagePointer
void SetInput2(const InputImageType *image)
Superclass::OutputImageRegionType OutputImageRegionType
const InputImageType * GetInput() const
InterpolateImageFunction< IntermediateImageType > InterpolatorType
Superclass::InputImageType InputImageType
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
Base class for all image interpolaters.
void SetInput1(const InputImageType *image)
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId) override
Superclass::InputImagePointer InputImagePointer
Base class for filters that take an image as input and produce an image as output.
TInputImage::PixelType InputPixelType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
void PrintSelf(std::ostream &os, Indent indent) const override
Superclass::OutputImageRegionType OutputImageRegionType
ImageToImageFilter< TInputImage, TOutputImage > Superclass
TOutputImage OutputImageType
InterpolatorType::Pointer m_Interpolator
#define itkConceptMacro(name, concept)
Interpolate an image from two N-D images.
IntermediateImageType::Pointer m_IntermediateImage
Templated n-dimensional image class.
Definition: itkImage.h:75
const InputImageType * GetInput2()