ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkBSplineResampleImageFilterBase.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 /*=========================================================================
19  *
20  * Portions of this file are subject to the VTK Toolkit Version 3 copyright.
21  *
22  * Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
23  *
24  * For complete copyright, license and disclaimer of warranty information
25  * please refer to the NOTICE file at the top of the ITK source tree.
26  *
27  *=========================================================================*/
28 #ifndef itkBSplineResampleImageFilterBase_h
29 #define itkBSplineResampleImageFilterBase_h
30 
31 #include <vector>
32 
34 #include "itkImageRegionIterator.h" // Used for the output iterator needs to
35  // match filter program
36 #include "itkProgressReporter.h"
37 #include "itkImageToImageFilter.h"
38 
39 namespace itk
40 {
81 template< typename TInputImage, typename TOutputImage >
82 class ITK_TEMPLATE_EXPORT BSplineResampleImageFilterBase:
83  public ImageToImageFilter< TInputImage, TOutputImage >
84 {
85 public:
86  ITK_DISALLOW_COPY_AND_ASSIGN(BSplineResampleImageFilterBase);
87 
93 
96 
98  // Must be sustantiated through another class. itkNewMacro( Self );
99 
101  using InputImageType = typename Superclass::InputImageType;
102 
104  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
105 
108 
110  using SizeType = typename TInputImage::SizeType;
111 
114 
116  using OutputImagePixelType = typename Superclass::OutputImagePixelType;
117 
120 
123 
126 
129  void SetSplineOrder(int SplineOrder);
130 
132  itkGetConstMacro(SplineOrder, int);
133 
134 protected:
136  void ReduceNDImage(OutputImageIterator & OutItr);
137 
139  void ExpandNDImage(OutputImageIterator & OutItr);
140 
143  virtual void InitializePyramidSplineFilter(int SplineOrder);
144 
146  virtual void Reduce1DImage(
147  const std::vector< double > & In,
148  OutputImageIterator & Iter,
149  unsigned int traverseSize,
150  ProgressReporter & progress
151  );
152 
154  virtual void Expand1DImage(
155  const std::vector< double > & In,
156  OutputImageIterator & Iter,
157  unsigned int traverseSize,
158  ProgressReporter & progress
159  );
160 
162  ~BSplineResampleImageFilterBase() override = default;
163  void PrintSelf(std::ostream & os, Indent indent) const override;
164 
165  int m_SplineOrder; // User specified spline order
166  int m_GSize; // downsampling filter size
167  int m_HSize; // upsampling filter size
168 
169  std::vector< double > m_G; // downsampling filter coefficients
170  std::vector< double > m_H; // upsampling filter coefficients
171 
172 private:
173 
174  // Resizes m_Scratch Variable based on image sizes
175  void InitializeScratch(SizeType DataLength);
176 
177  // Copies a line of data from the input to the m_Scratch for subsequent
178  // processing
179  void CopyInputLineToScratch(ConstInputImageIterator & Iter);
180 
181  void CopyOutputLineToScratch(ConstOutputImageIterator & Iter);
182 
183  void CopyLineToScratch(ConstInputImageIterator & Iter);
184 
185  std::vector< double > m_Scratch; // temp storage for processing
186  // of Coefficients
187 
188 };
189 } // namespace itk
190 
191 #ifndef ITK_MANUAL_INSTANTIATION
192 #include "itkBSplineResampleImageFilterBase.hxx"
193 #endif
194 
195 #endif
A multi-dimensional image iterator that visits image pixels within a region in a &quot;scan-line&quot; order...
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Base class for all process objects that output image data.
A multi-dimensional image iterator that visits image pixels within a region in a &quot;scan-line&quot; order...
typename OutputImageType::PixelType OutputImagePixelType
Uses the &quot;l2&quot; spline pyramid implementation of B-Spline Filters to up/down sample an image by a facto...
Implements progress tracking for a filter.
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49