ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkBSplineInterpolationWeightFunction.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 itkBSplineInterpolationWeightFunction_h
19#define itkBSplineInterpolationWeightFunction_h
20
21#include "itkFunctionBase.h"
22#include "itkContinuousIndex.h"
23#include "itkArray.h"
24#include "itkArray2D.h"
25#include "itkIndexRange.h"
26#include "itkMath.h"
27
28namespace itk
29{
47template <typename TCoordinate = float, unsigned int VSpaceDimension = 2, unsigned int VSplineOrder = 3>
48class ITK_TEMPLATE_EXPORT BSplineInterpolationWeightFunction
49 : public FunctionBase<ContinuousIndex<TCoordinate, VSpaceDimension>,
50 FixedArray<double, Math::UnsignedPower(VSplineOrder + 1, VSpaceDimension)>>
51{
52public:
53 ITK_DISALLOW_COPY_AND_MOVE(BSplineInterpolationWeightFunction);
54
58 FixedArray<double, Math::UnsignedPower(VSplineOrder + 1, VSpaceDimension)>>;
59
62
64 itkNewMacro(Self);
65
67 itkOverrideGetNameOfClassMacro(BSplineInterpolationWeightFunction);
68
70 static constexpr unsigned int SpaceDimension = VSpaceDimension;
71
73 static constexpr unsigned int SplineOrder = VSplineOrder;
74
77
79 static constexpr unsigned int NumberOfWeights{ WeightsType::Length };
80
84
87
89 static constexpr SizeType SupportSize{ SizeType::Filled(VSplineOrder + 1) };
90
93 WeightsType
94 Evaluate(const ContinuousIndexType & index) const override;
95
104 virtual void
105 Evaluate(const ContinuousIndexType & index, WeightsType & weights, IndexType & startIndex) const;
106
107#if !defined(ITK_LEGACY_REMOVE)
109 itkLegacyMacro(SizeType GetSupportSize() const)
110 {
111 return Self::SupportSize;
112 }
113
115 itkLegacyMacro(unsigned int GetNumberOfWeights() const)
116 {
117 return Self::NumberOfWeights;
118 }
119#endif
120
121protected:
124};
125
126} // end namespace itk
127
128#ifndef ITK_MANUAL_INSTANTIATION
129# include "itkBSplineInterpolationWeightFunction.hxx"
130#endif
131
132#endif
FunctionBase< ContinuousIndex< TCoordinate, VSpaceDimension >, FixedArray< double, Math::UnsignedPower(VSplineOrder+1, VSpaceDimension)> > Superclass
~BSplineInterpolationWeightFunction() override=default
virtual void Evaluate(const ContinuousIndexType &index, WeightsType &weights, IndexType &startIndex) const
ContinuousIndex< TCoordinate, VSpaceDimension > ContinuousIndexType
WeightsType Evaluate(const ContinuousIndexType &index) const override
A templated class holding a point in n-Dimensional image space.
Simulate a standard C array with copy semantics.
Implements transparent reference counting.
constexpr TReturnType UnsignedPower(const uintmax_t base, const uintmax_t exponent) noexcept
Definition itkMath.h:791
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Represent a n-dimensional index in a n-dimensional image.
Definition itkIndex.h:69
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition itkSize.h:70
static constexpr Self Filled(const SizeValueType value)
Definition itkSize.h:434