ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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{
70template <typename TInputImage, typename TOutputImage>
71class ITK_TEMPLATE_EXPORT BSplineResampleImageFilterBase : public ImageToImageFilter<TInputImage, TOutputImage>
72{
73public:
74 ITK_DISALLOW_COPY_AND_MOVE(BSplineResampleImageFilterBase);
75
81
83 itkOverrideGetNameOfClassMacro(BSplineResampleImageFilterBase);
84
86 // Must be instantiated through another class. itkNewMacro( Self );
87
89 using typename Superclass::InputImageType;
90
92 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
93
95 using IndexType = typename TInputImage::IndexType;
96
98 using SizeType = typename TInputImage::SizeType;
99
101 using RegionType = typename TInputImage::RegionType;
102
105
108
111
114
117 void
118 SetSplineOrder(int splineOrder);
119
121 itkGetConstMacro(SplineOrder, int);
122
123protected:
125 void
127
129 void
131
134 virtual void
136
142 virtual void
143 Reduce1DImage(const std::vector<double> & in,
145 unsigned int inTraverseSize,
146 ProgressReporter & progress);
147
153 virtual void
154 Expand1DImage(const std::vector<double> & in,
156 unsigned int inTraverseSize,
157 ProgressReporter & progress);
158
161 void
162 PrintSelf(std::ostream & os, Indent indent) const override;
163
164 int m_SplineOrder{}; // User specified spline order
165 int m_GSize{}; // downsampling filter size
166 int m_HSize{}; // upsampling filter size
167
168 std::vector<double> m_G{}; // downsampling filter coefficients
169 std::vector<double> m_H{}; // upsampling filter coefficients
170
171private:
173 void
175
177 void
179
180 void
182
183 void
185
186 std::vector<double> m_Scratch{}; // temp storage for processing
187 // of Coefficients
188};
189} // namespace itk
190
191#ifndef ITK_MANUAL_INSTANTIATION
192# include "itkBSplineResampleImageFilterBase.hxx"
193#endif
194
195#endif
void CopyInputLineToScratch(ConstInputImageIterator &Iter)
void CopyLineToScratch(ConstInputImageIterator &Iter)
itk::ImageLinearIteratorWithIndex< TOutputImage > OutputImageIterator
void ReduceNDImage(OutputImageIterator &outItr)
itk::ImageLinearConstIteratorWithIndex< TOutputImage > ConstOutputImageIterator
~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)
ImageToImageFilter< TInputImage, TOutputImage > Superclass
void PrintSelf(std::ostream &os, Indent indent) const override
void CopyOutputLineToScratch(ConstOutputImageIterator &Iter)
virtual void InitializePyramidSplineFilter(int SplineOrder)
void ExpandNDImage(OutputImageIterator &outItr)
itk::ImageLinearConstIteratorWithIndex< TInputImage > ConstInputImageIterator
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.
typename OutputImageType::PixelType OutputImagePixelType
Control indentation during Print() invocation.
Definition itkIndent.h:50
Implements progress tracking for a filter.
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....