ITK  6.0.0
Insight Toolkit
itkBSplineResampleImageFilterBase.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/*=========================================================================
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"
38
39namespace itk
40{
82template <typename TInputImage, typename TOutputImage>
83class ITK_TEMPLATE_EXPORT BSplineResampleImageFilterBase : public ImageToImageFilter<TInputImage, TOutputImage>
84{
85public:
86 ITK_DISALLOW_COPY_AND_MOVE(BSplineResampleImageFilterBase);
87
93
95 itkOverrideGetNameOfClassMacro(BSplineResampleImageFilterBase);
96
98 // Must be instantiated through another class. itkNewMacro( Self );
99
101 using typename Superclass::InputImageType;
102
104 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
105
108
111
114
116 using typename Superclass::OutputImagePixelType;
117
120
123
126
129 void
130 SetSplineOrder(int splineOrder);
131
133 itkGetConstMacro(SplineOrder, int);
134
135protected:
137 void
139
141 void
143
146 virtual void
148
154 virtual void
155 Reduce1DImage(const std::vector<double> & in,
157 unsigned int inTraverseSize,
158 ProgressReporter & progress);
159
165 virtual void
166 Expand1DImage(const std::vector<double> & in,
168 unsigned int inTraverseSize,
169 ProgressReporter & progress);
170
173 void
174 PrintSelf(std::ostream & os, Indent indent) const override;
175
176 int m_SplineOrder{}; // User specified spline order
177 int m_GSize{}; // downsampling filter size
178 int m_HSize{}; // upsampling filter size
179
180 std::vector<double> m_G{}; // downsampling filter coefficients
181 std::vector<double> m_H{}; // upsampling filter coefficients
182
183private:
185 void
187
189 void
191
192 void
194
195 void
197
198 std::vector<double> m_Scratch{}; // temp storage for processing
199 // of Coefficients
200};
201} // namespace itk
202
203#ifndef ITK_MANUAL_INSTANTIATION
204# include "itkBSplineResampleImageFilterBase.hxx"
205#endif
206
207#endif
Uses the "l2" spline pyramid implementation of B-Spline Filters to up/down sample an image by a facto...
void CopyInputLineToScratch(ConstInputImageIterator &Iter)
void CopyLineToScratch(ConstInputImageIterator &Iter)
void ReduceNDImage(OutputImageIterator &outItr)
~BSplineResampleImageFilterBase() override=default
void SetSplineOrder(int splineOrder)
virtual void Reduce1DImage(const std::vector< double > &in, OutputImageIterator &out, unsigned int inTraverseSize, ProgressReporter &progress)
void InitializeScratch(SizeType DataLength)
virtual void Expand1DImage(const std::vector< double > &in, OutputImageIterator &out, unsigned int inTraverseSize, ProgressReporter &progress)
void PrintSelf(std::ostream &os, Indent indent) const override
void CopyOutputLineToScratch(ConstOutputImageIterator &Iter)
virtual void InitializePyramidSplineFilter(int SplineOrder)
void ExpandNDImage(OutputImageIterator &outItr)
A multi-dimensional image iterator that visits image pixels within a region in a "scan-line" order.
A multi-dimensional image iterator that visits image pixels within a region in a "scan-line" order.
Base class for all process objects that output image data.
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Implements progress tracking for a filter.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....