Loading [MathJax]/jax/input/TeX/config.js
ITK 6.0.0
Insight Toolkit
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
Loading...
Searching...
No Matches
itkBSplineControlPointImageFilter.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#ifndef itkBSplineControlPointImageFilter_h
19#define itkBSplineControlPointImageFilter_h
20
24#include "itkFixedArray.h"
25#include "itkPointSet.h"
27#include "itkVector.h"
28#include "itkVectorContainer.h"
29
30#include "vnl/vnl_matrix.h"
31#include "vnl/vnl_vector.h"
32
33namespace itk
34{
35
55
56template <typename TInputImage, typename TOutputImage = TInputImage>
57class ITK_TEMPLATE_EXPORT BSplineControlPointImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
58{
59public:
60 ITK_DISALLOW_COPY_AND_MOVE(BSplineControlPointImageFilter);
61
66
68 itkNewMacro(Self);
69
71 itkOverrideGetNameOfClassMacro(BSplineControlPointImageFilter);
72
74 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
75
76 using ControlPointLatticeType = TInputImage;
77 using OutputImageType = TOutputImage;
78
80 using PixelType = typename OutputImageType::PixelType;
81 using RegionType = typename OutputImageType::RegionType;
82 using IndexType = typename OutputImageType::IndexType;
83 using PointType = typename OutputImageType::PointType;
84 using OutputImageRegionType = typename OutputImageType::RegionType;
85
86 using SpacingType = typename OutputImageType::SpacingType;
87 using OriginType = typename OutputImageType::PointType;
88 using SizeType = typename OutputImageType::SizeType;
89 using DirectionType = typename OutputImageType::DirectionType;
90
92 using RealType = float;
95
98
105
112
117 void
118 SetSplineOrder(unsigned int);
119
125
129 itkGetConstReferenceMacro(SplineOrder, ArrayType);
130
148 itkSetMacro(CloseDimension, ArrayType);
149 itkGetConstReferenceMacro(CloseDimension, ArrayType);
155 itkSetMacro(Spacing, SpacingType);
156 itkGetConstMacro(Spacing, SpacingType);
162 itkSetMacro(Origin, OriginType);
163 itkGetConstMacro(Origin, OriginType);
169 itkSetMacro(Size, SizeType);
170 itkGetConstMacro(Size, SizeType);
186 itkSetMacro(Direction, DirectionType);
187
191 itkGetConstMacro(Direction, DirectionType);
192
200 typename ControlPointLatticeType::Pointer RefineControlPointLattice(ArrayType);
201
202protected:
205 void
206 PrintSelf(std::ostream & os, Indent indent) const override;
207
209 void
211
212
213private:
218 void
220
225 unsigned int
226 SplitRequestedRegion(unsigned int, unsigned int, OutputImageRegionType &) override;
227
232 void
234
239
245
246 bool m_DoMultilevel{ false };
247 unsigned int m_MaximumNumberOfLevels{ 1 };
252
254
260
261 RealType m_BSplineEpsilon{ static_cast<RealType>(1e-3) };
262
263 inline typename RealImageType::IndexType
264 NumberToIndex(unsigned int number, typename RealImageType::SizeType size)
265 {
266 typename RealImageType::IndexType k;
267 k[0] = 1;
268
269 for (unsigned int i = 1; i < ImageDimension; ++i)
270 {
271 k[i] = size[ImageDimension - i - 1] * k[i - 1];
272 }
273 typename RealImageType::IndexType index;
274 for (unsigned int i = 0; i < ImageDimension; ++i)
275 {
276 index[ImageDimension - i - 1] = static_cast<unsigned int>(number / k[ImageDimension - i - 1]);
277 number %= k[ImageDimension - i - 1];
278 }
279 return index;
280 }
281};
282
283} // end namespace itk
284
285#ifndef ITK_MANUAL_INSTANTIATION
286# include "itkBSplineControlPointImageFilter.hxx"
287#endif
288
289#endif
FixedArray< unsigned int, Self::ImageDimension > ArrayType
void BeforeThreadedGenerateData() override
void PrintSelf(std::ostream &os, Indent indent) const override
vnl_matrix< RealType > m_RefinedLatticeCoefficients[ImageDimension]
typename OutputImageType::SpacingType SpacingType
unsigned int SplitRequestedRegion(unsigned int, unsigned int, OutputImageRegionType &) override
void CollapsePhiLattice(PointDataImageType *, PointDataImageType *, const RealType, const unsigned int)
typename OutputImageType::DirectionType DirectionType
typename OutputImageType::PointType PointType
~BSplineControlPointImageFilter() override=default
FixedArray< RealType, Self::ImageDimension > RealArrayType
RealImageType::IndexType NumberToIndex(unsigned int number, typename RealImageType::SizeType size)
typename OutputImageType::RegionType OutputImageRegionType
typename OutputImageType::PointType OriginType
Image< PointDataType, Self::ImageDimension > PointDataImageType
ControlPointLatticeType::Pointer RefineControlPointLattice(ArrayType)
typename PointSetType::PointDataContainer PointDataContainerType
typename OutputImageType::PixelType PixelType
void DynamicThreadedGenerateData(const OutputImageRegionType &) override
typename OutputImageType::RegionType RegionType
ImageToImageFilter< TInputImage, TOutputImage > Superclass
PointSet< PixelType, Self::ImageDimension > PointSetType
typename OutputImageType::IndexType IndexType
Image< RealType, Self::ImageDimension > RealImageType
typename PointDataImageType::Pointer PointDataImagePointer
BSpline kernel used for density estimation and nonparametric regression.
BSpline kernel used for density estimation and nonparametric regression.
Simulate a standard C array with copy semantics.
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition itkPointSet.h:82
typename MeshTraits::PixelType PixelType
typename MeshTraits::PointDataContainer PointDataContainer
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition itkSize.h:70